31#include <vtkSmartPointer.h> 
   32#include <vtkXMLReader.h> 
   33#include <vtkXMLWriter.h> 
   34#include <vtkXMLUnstructuredGridReader.h> 
   35#include <vtkXMLUnstructuredGridWriter.h> 
   36#include <vtkXMLPUnstructuredGridWriter.h> 
   37#include <vtkXMLPolyDataReader.h> 
   38#include <vtkXMLPolyDataWriter.h> 
   39#include <vtkXMLPPolyDataWriter.h> 
   40#include <vtkDataSetReader.h> 
   41#include <vtkDataSetWriter.h> 
   42#include <vtkPDataSetWriter.h> 
   43#include <vtkDataSet.h> 
   44#include <vtkUnstructuredGrid.h> 
   45#include <vtkPolyData.h> 
   46#include <vtkStructuredGrid.h> 
   47#include <vtkPointData.h> 
   48#include <vtkCellData.h> 
   49#include <vtkFieldData.h> 
   50#include <vtkCellTypes.h> 
   51#include <vtkDataArray.h> 
   52#include <vtkFloatArray.h> 
   53#include <vtkDoubleArray.h> 
   56#include <vtkTriangle.h> 
   60#include <vtkCellArray.h> 
   81  *ca = vtkSmartPointer< vtkCellArray >::New();
 
  122    w->SetDataModeToBinary();
 
  125    w->SetDataModeToAscii();
 
  136    w->SetDataModeToBinary();
 
  139    w->SetDataModeToAscii();
 
  169template <
class T, 
class TWriter, 
class PWriter>
 
  171                       const char *mfilename,
int metricData,
 
  172                       int binary,
int npart,
int myid,
int master) {
 
  173  int hasPointRef = 0, hasCellRef = 0;
 
  176  vtkSmartPointer<vtkPoints> points = vtkSmartPointer<vtkPoints>::New();
 
  179  for ( MMG5_int k=1; k<=
mesh->
np; ++k ) {
 
  181    if ( !
MG_VOK(ppt) ) 
continue;
 
  185    if ( 
mesh->
dim == 2 ) ppt->
c[2] = 0.;
 
  186    if ( ppt->
ref )       hasPointRef = 1;
 
  188    points->InsertNextPoint(ppt->
c[0], ppt->
c[1], ppt->
c[2]);
 
  192  auto dataset = vtkSmartPointer<T> :: New();
 
  193  dataset->SetPoints(points);
 
  197  vtkSmartPointer<vtkCellArray>  cellArray  = vtkSmartPointer<vtkCellArray>  ::New();
 
  198  vtkSmartPointer<vtkWedge>      wedge      = vtkSmartPointer<vtkWedge>      ::New();
 
  199  vtkSmartPointer<vtkTetra>      tetra      = vtkSmartPointer<vtkTetra>      ::New();
 
  200  vtkSmartPointer<vtkTriangle>   tria       = vtkSmartPointer<vtkTriangle>   ::New();
 
  201  vtkSmartPointer<vtkQuad>       quadra     = vtkSmartPointer<vtkQuad>       ::New();
 
  202  vtkSmartPointer<vtkLine>       edge       = vtkSmartPointer<vtkLine>       ::New();
 
  211  for ( MMG5_int k=1; k<=
mesh->
na; ++k ) {
 
  213    if ( !pa || !pa->
a ) 
continue;
 
  215    if ( pa->
ref ) hasCellRef = 1;
 
  220    cellArray->InsertNextCell(edge);
 
  221    types[ic++] = VTK_LINE;
 
  227 for ( MMG5_int k=1; k<=
mesh->
nt; ++k ) {
 
  229   if ( !
MG_EOK(ptt) ) 
continue;
 
  231   if ( ptt->
ref ) hasCellRef = 1;
 
  233   for ( 
int i=0; i<3; ++i ) {
 
  236   cellArray->InsertNextCell(tria);
 
  237   types[ic++] = VTK_TRIANGLE;
 
  241 for ( MMG5_int k=1; k<=
mesh->
nquad; ++k ) {
 
  243   if ( !
MG_EOK(pq) ) 
continue;
 
  245   if ( pq->
ref ) hasCellRef = 1;
 
  247   for ( 
int i=0; i<4; ++i ) {
 
  250   cellArray->InsertNextCell(quadra);
 
  251   types[ic++] = VTK_QUAD;
 
  255  for ( MMG5_int k=1; k<=
mesh->
ne; ++k ) {
 
  257    if ( !
MG_EOK(pt) ) 
continue;
 
  259    if ( pt->
ref ) hasCellRef = 1;
 
  261    for ( 
int i=0; i<4; ++i ) {
 
  264    cellArray->InsertNextCell(tetra);
 
  265    types[ic++] = VTK_TETRA;
 
  271    if ( !
MG_EOK(ppr) ) 
continue;
 
  273    if ( ppr->
ref ) hasCellRef = 1;
 
  275    for ( 
int i=0; i<6; ++i ) {
 
  278    cellArray->InsertNextCell(wedge);
 
  279    types[ic++] = VTK_WEDGE;
 
  289    auto *ar = vtkIntArray::New();
 
  291    ar->SetNumberOfComponents(1);
 
  292    ar->SetNumberOfTuples(
mesh->
np);
 
  293    ar->SetName(
"medit:ref");
 
  295    for ( MMG5_int k = 0; k < 
mesh->
np; k++ ) {
 
  297      if ( !
MG_VOK(ppt) ) 
continue;
 
  299      ar->SetTuple1(ppt->
tmp,ppt->
ref);
 
  302    dataset->GetPointData()->AddArray(ar);
 
  305    auto *ar = vtkIntArray::New();
 
  307    ar->SetNumberOfComponents(1);
 
  308    ar->SetNumberOfTuples(nc);
 
  309    ar->SetName(
"medit:ref");
 
  312    for ( MMG5_int k = 1; k <= 
mesh->
na; k++ ) {
 
  316    for ( MMG5_int k = 1; k <= 
mesh->
nt; k++ ) {
 
  320    for ( MMG5_int k = 1; k <= 
mesh->
nquad; k++ ) {
 
  324    for ( MMG5_int k = 1; k <= 
mesh->
ne; k++ ) {
 
  328    for ( MMG5_int k = 1; k <= 
mesh->
nprism; k++ ) {
 
  333    dataset->GetCellData()->AddArray(ar);
 
  343  if ( metricData==1 ) {
 
  356  static int mmgWarn = 0;
 
  357  for ( 
int isol=0; isol<nsols; ++isol) {
 
  359    if (metricData && isol==0) {
 
  363    else if (metricData){
 
  364      psl = &
sol[1][isol-1];
 
  374        fprintf(stderr, 
"  ## Warning: %s: missing data for at least 1 solution." 
  375                " Skipped.\n",__func__);
 
  381    if ( psl->
size == 1 ) {
 
  384    else if ( psl->
size == psl->
dim ) {
 
  391    auto *ar = vtkDoubleArray::New();
 
  393    ar->SetNumberOfComponents(ncp);
 
  394    ar->SetNumberOfTuples(
mesh->
np);
 
  408      if ( metricData && strstr(data,
":ls")) {
 
  409        memmove(data2,data,strlen(data)-3);
 
  410        strcat ( data2 , 
":metric");
 
  413      else if (metricData && !strstr(data,
":metric") && isol==0) {
 
  414        strcat ( data , 
":metric");
 
  424      ar->SetName(
"no_name");
 
  430    if ( psl->
size!= (psl->
dim*(psl->
dim+1))/2 ) {
 
  432      for ( MMG5_int k=1; k<=
mesh->
np; k++) {
 
  434        if ( !
MG_VOK(ppt) ) 
continue;
 
  436        MMG5_int    iadr  = k*psl->
size;
 
  437        for ( 
int i=0; i<psl->
size; ++i ) {
 
  438          dfmt[i] = psl->
m[iadr+i];
 
  441        if ( psl->
dim==2 && ncp==3 ) {
 
  444        ar->SetTuple(ppt->
tmp,dfmt);
 
  449      for ( MMG5_int k=1; k<=
mesh->
np; k++) {
 
  451        if ( !
MG_VOK(ppt) ) 
continue;
 
  453        MMG5_int iadr  = k*psl->
size;
 
  454        double *d = &psl->
m[iadr];
 
  462          double dbuf[6] = {0,0,0,0,0,0};
 
  468            for ( 
int i=0; i<psl->
size; i++)  dbuf[i] = psl->
m[psl->
size*k+i];
 
  480        ar->SetTuple(ppt->
tmp,dfmt);
 
  483    dataset->GetPointData()->AddArray(ar);
 
  490    vtkSmartPointer<PWriter> writer = vtkSmartPointer<PWriter>::New();
 
  492#if VTK_MAJOR_VERSION <= 5 
  493    writer->SetInput(dataset);
 
  495    writer->SetInputData(dataset);
 
  498    writer->SetFileName(mfilename);
 
  500    writer->SetNumberOfPieces(npart);
 
  501    writer->SetStartPiece(myid);
 
  502    writer->SetEndPiece(myid);
 
  507    vtkSmartPointer<TWriter> writer = vtkSmartPointer<TWriter>::New();
 
  509    writer->SetFileName(mfilename);
 
  511#if VTK_MAJOR_VERSION <= 5 
  512    writer->SetInput(dataset);
 
  514    writer->SetInputData(dataset);
 
  538template <
class T, 
class TWriter,
class PWriter>
 
  540                     int metricData,
int binary) {
 
  543                                                metricData,binary,0,0,0);
 
char * MMG5_Get_basename(char *path)
 
MMG5_pMesh MMG5_pSol * sol
 
MMG5_pMesh char * filename
 
void MMG5_build3DMetric(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int ip, double dbuf[6])
 
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
 
#define MMG5_SAFE_MALLOC(ptr, size, type, law)
 
#define MMG5_SAFE_FREE(ptr)
 
#define MMG5_DEL_MEM(mesh, ptr)
 
Structure to store edges of am MMG mesh.
 
Structure to store vertices of an MMG mesh.
 
Structure to store prsim of a MMG mesh.
 
Structure to store quadrangles of an MMG mesh.
 
Structure to store tetrahedra of an MMG mesh.
 
Structure to store triangles of a MMG mesh.
 
int MMG5_loadVtuMesh_part1(MMG5_pMesh, const char *, vtkDataSet **, int8_t *, int8_t *, int *, int8_t *, int8_t *)
 
static void MMG5_internal_VTKSetLine(vtkSmartPointer< vtkPolyData > d, vtkSmartPointer< vtkCellArray > *ca)
 
int MMG5_loadVtpMesh_part1(MMG5_pMesh, const char *, vtkDataSet **, int8_t *, int8_t *, int *, int8_t *, int8_t *)
 
int MMG5_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename, int metricData, int binary)
 
static void MMG5_internal_VTKbinary(vtkXMLUnstructuredGridWriter *w, int binary)
 
int MMG5_saveVtkMesh_i(MMG5_pMesh mesh, MMG5_pSol *sol, const char *mfilename, int metricData, int binary, int npart, int myid, int master)
 
int MMG5_loadVtkMesh_part2(MMG5_pMesh, MMG5_pSol *, vtkDataSet **, int8_t, int8_t, int, int8_t, int8_t)
 
static void MMG5_internal_VTKSetCells(int *t, vtkSmartPointer< vtkPolyData > d, vtkSmartPointer< vtkCellArray > ca)
 
int MMG5_loadVtkMesh_part1(MMG5_pMesh, const char *, vtkDataSet **, int8_t *, int8_t *, int *, int8_t *, int8_t *)