32#include <vtkSmartPointer.h>
33#include <vtkXMLReader.h>
34#include <vtkXMLWriter.h>
35#include <vtkXMLUnstructuredGridReader.h>
36#include <vtkXMLPolyDataReader.h>
37#include <vtkDataSetReader.h>
38#include <vtkDataSet.h>
39#include <vtkUnstructuredGrid.h>
40#include <vtkPolyData.h>
41#include <vtkStructuredGrid.h>
42#include <vtkPointData.h>
43#include <vtkCellData.h>
44#include <vtkFieldData.h>
45#include <vtkCellTypes.h>
46#include <vtkDataArray.h>
62template<
class TReader>
65 vtkSmartPointer<TReader> reader =
66 vtkSmartPointer<TReader>::New();
67 reader->SetFileName(fileName);
70 if ( !reader->GetOutput() ) {
71 throw "Unable to open file.";
74 reader->GetOutput()->Register(reader);
76 return vtkDataSet::SafeDownCast(reader->GetOutput());
104 int8_t *ptMeditRef, int8_t *eltMeditRef,
105 int *nsols, int8_t *metricData, int8_t *lsData ) {
107 static int8_t mmgWarn0 = 0;
108 static int8_t mmgWarn1 = 0;
117 mesh->
np = dataset->GetNumberOfPoints();
118 vtkIdType numCells = dataset->GetNumberOfCells();
121 for ( vtkIdType i = 0; i < numCells; i++ ) {
122 int typ = dataset->GetCellType(i);
127 case ( VTK_POLY_LINE ):
128 printf(
"polylin %lld \n", dataset->GetCell(i)->GetNumberOfPoints() - 1 );
130 printf(
"%lld %lld %lld\n", dataset->GetCell(i)->GetPointId(0)+1,
131 dataset->GetCell(i)->GetPointId(1)+1,dataset->GetCell(i)->GetPointId(2)+1);
132 mesh->
na += dataset->GetCell(i)->GetNumberOfPoints() - 1;
137 case ( VTK_TRIANGLE ):
151 printf(
" ## Warning:%s: unexpected element type (%d).",__func__,typ);
160 auto *pd = dataset->GetPointData();
163 int nmetricField = 0;
168 npointData = pd->GetNumberOfArrays();
169 for (
int k = 0; k < npointData; k++) {
170 if ( strstr(pd->GetArrayName(k),
":metric" ) ) {
173 else if ( strstr(pd->GetArrayName(k),
":ls" ) ) {
176 else if ( strstr(pd->GetArrayName(k),
"medit:ref" ) ) {
182 if ( npointRef > 1 ) {
183 printf(
" ## Warning:%s: %d reference fields detected (labelled 'medit:ref')."
184 " Only the last is used, others are ignored.", __func__, npointRef);
186 if ( nmetricField > 1 ) {
187 printf(
" ## Error:%s: %d metric fields detected (labelled with a string"
188 " containing the 'metric' keyword).\n"
189 " Exit Program.\n",__func__,nmetricField);
192 if ( nlsField > 1 ) {
193 printf(
" ## Error:%s: %d ls fields detected (labelled with a string"
194 " containing the 'ls' keyword).\n"
195 " Exit Program.\n",__func__,nlsField);
201 auto *cd = dataset->GetCellData();
207 ncellData = cd->GetNumberOfArrays();
208 for (
int k = 0; k < ncellData; k++) {
209 if ( strstr(cd->GetArrayName(k),
"medit:ref" ) ) {
217 auto *fd = dataset->GetFieldData();
218 if ( fd->GetNumberOfArrays() ) {
219 printf(
" ## Warning:%s: VTK field data not used in Mmg."
220 " Ignored.\n",__func__ );
223 *nsols = npointData + ncellData - npointRef - ncellRef;
224 *metricData = ( nmetricField > 0 );
225 *lsData = ( nlsField > 0 );
246 int8_t *ptMeditRef,int8_t *eltMeditRef,
int *nsols,
247 int8_t *metricData, int8_t *lsData) {
250 (*ptMeditRef) = (*eltMeditRef) = -1;
254 (*dataset) = MMG5_load_vtkXMLFile<vtkXMLPolyDataReader> (
filename );
262 nsols,metricData,lsData);
286 int8_t *ptMeditRef,int8_t *eltMeditRef,
int *nsols,
287 int8_t *metricData, int8_t *lsData) {
290 (*ptMeditRef) = (*eltMeditRef) = -1;
294 (*dataset) = MMG5_load_vtkXMLFile<vtkDataSetReader> (
filename );
302 nsols,metricData,lsData);
326 int8_t *ptMeditRef,int8_t *eltMeditRef,
int *nsols,
327 int8_t *metricData, int8_t *lsData) {
330 (*ptMeditRef) = (*eltMeditRef) = -1;
335 (*dataset) = MMG5_load_vtkXMLFile<vtkXMLUnstructuredGridReader> (
filename );
343 nsols,metricData,lsData);
367 int8_t ptMeditRef,int8_t eltMeditRef,
int nsols,
368 int8_t metricData, int8_t lsData) {
369 vtkSmartPointer<vtkDataArray> ptar = NULL, car = NULL;
372 static int8_t mmgWarn1 = 0;
375 if ( ptMeditRef > -1 ) {
376 auto *pd = (*dataset)->GetPointData();
377 ptar = pd->GetArray(ptMeditRef);
380 assert ( strstr( ptar->GetName(),
"medit:ref" ) );
383 assert ( ptar->GetNumberOfComponents() == 1 );
385 MMG5_int np = ptar->GetNumberOfTuples();
387 printf(
" ## Error: Point data size (%" MMG5_PRId
") differs from the number of"
388 " vertices (%" MMG5_PRId
")\n",np,
mesh->
np);
392 for ( vtkIdType k = 0; k < (*dataset)->GetNumberOfPoints(); k++ ) {
397 ppt->
ref = ptar->GetTuple1(k);
398 if ( ppt->
ref < 0 ) {
406 for ( vtkIdType k = 0; k < (*dataset)->GetNumberOfPoints(); k++ ) {
431 vtkIdType numCells = (*dataset)->GetNumberOfCells();
433 if ( eltMeditRef > -1 ) {
435 car = (*dataset)->GetCellData()->GetArray(eltMeditRef);
438 assert ( strstr( car->GetName(),
"medit:ref" ) );
440 assert ( car->GetNumberOfComponents() == 1 );
442 ne = car->GetNumberOfTuples();
443 if ( ne != numCells ) {
444 printf(
" ## Error: Cell data size (%" MMG5_PRId
") differs from the number of"
445 " cells (%lld)\n",ne,numCells);
451 for ( vtkIdType k = 0; k < numCells; k++ ) {
459 int typ = (*dataset)->GetCellType(k);
465 ppt->
ref = car ? car->GetTuple1(k) : 0;
466 if ( ppt->
ref < 0 ) {
472 case ( VTK_POLY_LINE ):
474 n = (*dataset)->GetCell(k)->GetNumberOfPoints() - 1;
476 for (
int i=0; i<n; ++i ) {
478 ref = car ? car->GetTuple1(k) : 0;
486 for (
int i=0; i < n; ++i ) {
488 pa->
a = (*dataset)->GetCell(k)->GetPointId(i) +1;
489 pa->
b = (*dataset)->GetCell(k)->GetPointId(i+1)+1;
502 ref = car ? car->GetTuple1(k) : 0;
511 pa->
a = (*dataset)->GetCell(k)->GetPointId(0)+1;
512 pa->
b = (*dataset)->GetCell(k)->GetPointId(1)+1;
520 assert( na+nbl_a <= mesh->na );
523 case ( VTK_TRIANGLE ):
525 ref = car ? car->GetTuple1(k) : 0;
534 for (
int i=0; i<3; ++i ) {
535 ptt->
v[i] = (*dataset)->GetCell(k)->GetPointId(i)+1;
539 if ( ptt->
ref < 0 ) {
548 pq->
ref = car ? car->GetTuple1(k) : 0;
549 for (
int i=0; i<4; ++i ) {
550 pq->
v[i] = (*dataset)->GetCell(k)->GetPointId(i)+1;
563 pt->
ref = car ? car->GetTuple1(k) : 0;
564 for (
int i=0; i<4; ++i ) {
565 pt->
v[i] = (*dataset)->GetCell(k)->GetPointId(i)+1;
578 ppr->
ref = car ? car->GetTuple1(k) : 0;
579 for (
int i=0; i<6; ++i ) {
580 ppr->
v[i] = (*dataset)->GetCell(k)->GetPointId(i)+1;
582 if ( ppr->
ref < 0 ) {
590 printf(
" ## Warning:%s: unexpected element type (%d).\n",
603 assert ( na + nbl_a ==
mesh->
na );
604 assert ( nt + nbl_t ==
mesh->
nt );
612 else if ( nt < mesh->nt ) {
614 fprintf(stderr,
" Exit program.\n");
624 else if ( na < mesh->na ) {
626 fprintf(stderr,
" Exit program.\n");
636 if (
ier < 1 )
return ier;
649 auto *pd = (*dataset)->GetPointData();
650 auto *cd = (*dataset)->GetCellData();
653 int npointData = pd->GetNumberOfArrays();
655 for (
int j = 0; j < npointData; j++) {
656 bool metricField = 0;
658 strcpy(chaine,pd->GetArrayName(j));
660 if ( strstr(chaine,
"medit:ref" ) ) {
663 else if ( strstr(chaine,
":metric") ) {
667 if ( metricField && (metricData*lsData) ) {
683 fprintf(stderr,
"\n ## Warning: %s: unable to set solution name for"
684 " at least 1 solution.\n",__func__);
692 fprintf(stderr,
"\n ## Warning: %s: unable to set solution name for"
693 " at least 1 solution.\n",__func__);
698 auto ar = pd->GetArray(j);
700 psl->
np = ar->GetNumberOfTuples();
702 fprintf(stderr,
" ** MISMATCHES DATA: THE NUMBER OF VERTICES IN "
703 "THE MESH (%" MMG5_PRId
") DIFFERS FROM THE NUMBER OF VERTICES IN "
704 "THE SOLUTION (%" MMG5_PRId
") \n",
mesh->
np,psl->
np);
708 int ncp = ar->GetNumberOfComponents();
714 else if ( ncp == 2 || ncp == 3 ) {
730 fprintf(stderr,
" ** UNEXPECTED NUMBER OF COMPONENTS (%d). IGNORED \n",ncp);
739 fprintf(stderr,
" Exit program.\n");
743 switch ( psl->
type ) {
744 case ( 1 ):
case ( 2 ):
745 for (MMG5_int k=1; k<=psl->
np; k++) {
746 MMG5_int iadr = k*psl->
size;
747 ar->GetTuple ( k-1, &psl->
m[iadr] );
755 for (MMG5_int k=1; k<=psl->
np; k++) {
756 ar->GetTuple ( k-1, dbuf );
757 MMG5_int iadr = psl->
size*k;
761 if ( psl->
dim ==2 ) {
762 psl->
m[iadr] = dbuf[0];
763 psl->
m[iadr+1] = dbuf[1];
764 psl->
m[iadr+2] = dbuf[3];
765 psl->
m[iadr+3] = dbuf[4];
768 for (
int i=0 ; i<9 ; i++ ) {
769 psl->
m[iadr+i] = dbuf[i];
775 if ( psl->
dim ==2 ) {
776 assert ( dbuf[1] == dbuf[2] );
778 psl->
m[iadr] = dbuf[0];
779 psl->
m[iadr+1] = dbuf[1];
780 psl->
m[iadr+2] = dbuf[3];
783 assert ( dbuf[1]==dbuf[3] || dbuf[2]==dbuf[6] || dbuf[5]==dbuf[7] );
785 psl->
m[iadr+0] = dbuf[0];
786 psl->
m[iadr+1] = dbuf[1];
787 psl->
m[iadr+2] = dbuf[2];
788 psl->
m[iadr+3] = dbuf[4];
789 psl->
m[iadr+4] = dbuf[5];
790 psl->
m[iadr+5] = dbuf[8];
797 fprintf(stderr,
" ** UNEXPECTED METRIC TYPE (%d). EXIT PROGRAM \n",psl->
type);
804 int ncellData = cd->GetNumberOfArrays();
806 for (
int j = 0; j < ncellData; j++) {
809 strcpy(chaine,cd->GetArrayName(j));
811 if ( strstr(chaine,
"medit:ref" ) ) {
824 fprintf(stderr,
"\n ## Warning: %s: unable to set solution name for"
825 " at least 1 solution.\n",__func__);
829 auto ar = cd->GetArray(j);
831 psl->
np = ar->GetNumberOfTuples();
832 if ( numCells != psl->
np ) {
833 fprintf(stderr,
" ** MISMATCHES DATA: THE NUMBER OF ELEMENTS IN "
834 "THE MESH (%lld) DIFFERS FROM THE NUMBER OF CELLS DATA IN "
835 "THE SOLUTION (%" MMG5_PRId
") \n",numCells,psl->
np);
839 int ncp = ar->GetNumberOfComponents();
845 else if ( ncp == 2 || ncp == 3 ) {
855 fprintf(stderr,
" ** UNEXPECTED NUMBER OF COMPONENTS (%d). IGNORED \n",ncp);
864 fprintf(stderr,
" Exit program.\n");
868 switch ( psl->
type ) {
869 case ( 1 ):
case ( 2 ):
870 for (MMG5_int k=1; k<=psl->
np; k++) {
871 MMG5_int iadr = k*psl->
size;
872 ar->GetTuple ( k-1, &psl->
m[iadr] );
880 for (MMG5_int k=1; k<=psl->
np; k++) {
881 ar->GetTuple ( k-1, dbuf );
882 MMG5_int iadr = psl->
size*k;
885 if ( psl->
dim ==2 ) {
886 psl->
m[iadr] = dbuf[0];
887 psl->
m[iadr+1] = dbuf[1];
888 psl->
m[iadr+2] = dbuf[3];
889 psl->
m[iadr+3] = dbuf[4];
892 for (
int i=0 ; i<9 ; i++ ) {
893 psl->
m[iadr+i] = dbuf[i];
900 fprintf(stderr,
" ** UNEXPECTED METRIC TYPE (%d). EXIT PROGRAM \n",psl->
type);
911 (*dataset)->Delete();
int MMG5_Set_inputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solin)
MMG5_pMesh MMG5_pSol * sol
MMG5_pMesh char * filename
int MMG5_check_readedMesh(MMG5_pMesh mesh, MMG5_int nref)
Types used throughout the Mmg libraries.
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
#define MMG5_ADD_MEM(mesh, size, message, law)
#define MMG5_FILESTR_LGTH
#define MMG5_DEL_MEM(mesh, ptr)
#define MMG5_SAFE_RECALLOC(ptr, prevSize, newSize, type, message, law)
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_loadVtkMesh_part1(MMG5_pMesh mesh, const char *filename, vtkDataSet **dataset, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData, int8_t *lsData)
int MMG5_loadVtpMesh_part1(MMG5_pMesh mesh, const char *filename, vtkDataSet **dataset, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData, int8_t *lsData)
vtkDataSet * MMG5_load_vtkXMLFile(const char *fileName)
I/O at VTK format.
int MMG5_loadVtuMesh_part1(MMG5_pMesh mesh, const char *filename, vtkDataSet **dataset, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData, int8_t *lsData)
static int MMG5_count_vtkEntities(vtkDataSet *dataset, MMG5_pMesh mesh, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData, int8_t *lsData)
int MMG5_loadVtkMesh_part2(MMG5_pMesh mesh, MMG5_pSol *sol, vtkDataSet **dataset, int8_t ptMeditRef, int8_t eltMeditRef, int nsols, int8_t metricData, int8_t lsData)