Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
vtkparser.cpp
Go to the documentation of this file.
1/* =============================================================================
2** This file is part of the mmg software package for the tetrahedral
3** mesh modification.
4** Copyright (c) Bx INP/CNRS/Inria/UBordeaux/UPMC, 2004-
5**
6** mmg is free software: you can redistribute it and/or modify it
7** under the terms of the GNU Lesser General Public License as published
8** by the Free Software Foundation, either version 3 of the License, or
9** (at your option) any later version.
10**
11** mmg is distributed in the hope that it will be useful, but WITHOUT
12** ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13** FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
14** License for more details.
15**
16** You should have received a copy of the GNU Lesser General Public
17** License and of the GNU General Public License along with mmg (in
18** files COPYING.LESSER and COPYING). If not, see
19** <http://www.gnu.org/licenses/>. Please read their terms carefully and
20** use this copy of the mmg distribution only if you accept them.
21** =============================================================================
22*/
23
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>
47
48#include "vtkparser.hpp"
49#include "libmmgtypes.h"
50#include "mmgcommon_private.h"
51
62template<class TReader>
63vtkDataSet *MMG5_load_vtkXMLFile(const char*fileName)
64{
65 vtkSmartPointer<TReader> reader =
66 vtkSmartPointer<TReader>::New();
67 reader->SetFileName(fileName);
68
69 reader->Update();
70 if ( !reader->GetOutput() ) {
71 throw "Unable to open file.";
72 }
73
74 reader->GetOutput()->Register(reader);
75
76 return vtkDataSet::SafeDownCast(reader->GetOutput());
77}
78
102static
103int MMG5_count_vtkEntities ( vtkDataSet *dataset, MMG5_pMesh mesh,
104 int8_t *ptMeditRef, int8_t *eltMeditRef,
105 int *nsols, int8_t *metricData, int8_t *lsData ) {
106
107 static int8_t mmgWarn0 = 0;
108 static int8_t mmgWarn1 = 0;
109
110 // Count the number of entities in the mesh
111 assert ( !mesh->np );
112 assert ( !mesh->nt );
113 assert ( !mesh->na );
114 assert ( !mesh->nquad );
115 assert ( !mesh->nprism );
116
117 mesh->np = dataset->GetNumberOfPoints();
118 vtkIdType numCells = dataset->GetNumberOfCells();
119
120 assert ( !mesh->npi );
121 for ( vtkIdType i = 0; i < numCells; i++ ) {
122 int typ = dataset->GetCellType(i);
123 switch ( typ ) {
124 case ( VTK_VERTEX ):
125 ++mesh->npi;
126 break;
127 case ( VTK_POLY_LINE ):
128 printf( "polylin %lld \n", dataset->GetCell(i)->GetNumberOfPoints() - 1 );
129
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;
133 break;
134 case ( VTK_LINE ):
135 ++mesh->na;
136 break;
137 case ( VTK_TRIANGLE ):
138 ++mesh->nt;
139 break;
140 case ( VTK_QUAD):
141 ++mesh->nquad;
142 break;
143 case ( VTK_TETRA):
144 ++mesh->ne;
145 break;
146 case ( VTK_WEDGE):
147 ++mesh->nprism;
148 break;
149 default:
150 if ( !mmgWarn1 ) {
151 printf( " ## Warning:%s: unexpected element type (%d).",__func__,typ);
152 mmgWarn1 = 1;
153 }
154 }
155 }
156 assert ( mesh->npi == mesh->np || !mesh->npi );
157 mesh->npi = 0;
158
159 // Count the number of point data and detect point references
160 auto *pd = dataset->GetPointData();
161 int npointData = 0;
162 int npointRef = 0;
163 int nmetricField = 0;
164 int nlsField = 0;
165
166 *ptMeditRef = -1;
167 if ( pd ) {
168 npointData = pd->GetNumberOfArrays();
169 for (int k = 0; k < npointData; k++) {
170 if ( strstr(pd->GetArrayName(k),":metric" ) ) {
171 ++nmetricField;
172 }
173 else if ( strstr(pd->GetArrayName(k),":ls" ) ) {
174 ++nlsField;
175 }
176 else if ( strstr(pd->GetArrayName(k),"medit:ref" ) ) {
177 (*ptMeditRef) = k;
178 ++npointRef;
179 }
180 }
181
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);
185 }
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);
190 return -1;
191 }
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);
196 return -1;
197 }
198 }
199
200 // Count the number of cell data and detect cell references
201 auto *cd = dataset->GetCellData();
202 int ncellData = 0;
203 int ncellRef = 0;
204
205 *eltMeditRef = -1;
206 if ( cd ) {
207 ncellData = cd->GetNumberOfArrays();
208 for (int k = 0; k < ncellData; k++) {
209 if ( strstr(cd->GetArrayName(k),"medit:ref" ) ) {
210 *eltMeditRef = k;
211 ++ncellRef;
212 }
213 }
214 }
215
216 // Count the number of field data
217 auto *fd = dataset->GetFieldData();
218 if ( fd->GetNumberOfArrays() ) {
219 printf( " ## Warning:%s: VTK field data not used in Mmg."
220 " Ignored.\n",__func__ );
221 }
222
223 *nsols = npointData + ncellData - npointRef - ncellRef;
224 *metricData = ( nmetricField > 0 );
225 *lsData = ( nlsField > 0 );
226
227 return 1;
228}
229
245int MMG5_loadVtpMesh_part1(MMG5_pMesh mesh,const char *filename,vtkDataSet **dataset,
246 int8_t *ptMeditRef,int8_t *eltMeditRef,int *nsols,
247 int8_t *metricData, int8_t *lsData) {
248
249 (*nsols) = 0;
250 (*ptMeditRef) = (*eltMeditRef) = -1;
251
252 // Read all the data from the file
253 try {
254 (*dataset) = MMG5_load_vtkXMLFile<vtkXMLPolyDataReader> ( filename );
255 }
256 catch( ... ) {
257 return 0;
258 }
259
260 // count the number of entities of each type
261 int ier = MMG5_count_vtkEntities ( (*dataset),mesh,ptMeditRef,eltMeditRef,
262 nsols,metricData,lsData);
263
264 if ( ier != 1 ) {
265 return -1;
266 }
267 return 1;
268}
269
285int MMG5_loadVtkMesh_part1(MMG5_pMesh mesh,const char *filename,vtkDataSet **dataset,
286 int8_t *ptMeditRef,int8_t *eltMeditRef,int *nsols,
287 int8_t *metricData, int8_t *lsData) {
288
289 (*nsols) = 0;
290 (*ptMeditRef) = (*eltMeditRef) = -1;
291
292 // Read all the data from the file
293 try {
294 (*dataset) = MMG5_load_vtkXMLFile<vtkDataSetReader> ( filename );
295 }
296 catch( ... ) {
297 return 0;
298 }
299
300 // count the number of entities of each type
301 int ier = MMG5_count_vtkEntities ( (*dataset),mesh,ptMeditRef,eltMeditRef,
302 nsols,metricData,lsData);
303
304 if ( ier != 1 ) {
305 return -1;
306 }
307 return 1;
308}
309
325int MMG5_loadVtuMesh_part1(MMG5_pMesh mesh,const char *filename,vtkDataSet **dataset,
326 int8_t *ptMeditRef,int8_t *eltMeditRef,int *nsols,
327 int8_t *metricData, int8_t *lsData) {
328
329 (*nsols) = 0;
330 (*ptMeditRef) = (*eltMeditRef) = -1;
331 (*metricData) = 0;
332
333 // Read all the data from the file
334 try {
335 (*dataset) = MMG5_load_vtkXMLFile<vtkXMLUnstructuredGridReader> ( filename );
336 }
337 catch ( ... ) {
338 return 0;
339 }
340
341 // count the number of entities of each type
342 int ier = MMG5_count_vtkEntities ( (*dataset),mesh,ptMeditRef,eltMeditRef,
343 nsols,metricData,lsData);
344
345 if ( ier != 1 ) {
346 return -1;
347 }
348
349 return 1;
350}
351
352
353
367 int8_t ptMeditRef,int8_t eltMeditRef,int nsols,
368 int8_t metricData, int8_t lsData) {
369 vtkSmartPointer<vtkDataArray> ptar = NULL, car = NULL;
370 int ier;
371 MMG5_int nref = 0;
372 static int8_t mmgWarn1 = 0;
373
374 // Point transfers in Mmg data structure
375 if ( ptMeditRef > -1 ) {
376 auto *pd = (*dataset)->GetPointData();
377 ptar = pd->GetArray(ptMeditRef);
378
379 // Check the field name
380 assert ( strstr( ptar->GetName(), "medit:ref" ) );
381
382 // Check that we get 1 data only
383 assert ( ptar->GetNumberOfComponents() == 1 );
384
385 MMG5_int np = ptar->GetNumberOfTuples();
386 if ( np != mesh->np ) {
387 printf( " ## Error: Point data size (%" MMG5_PRId ") differs from the number of"
388 " vertices (%" MMG5_PRId ")\n",np,mesh->np);
389 return -1;
390 }
391 // read vertices and vertices refs
392 for ( vtkIdType k = 0; k < (*dataset)->GetNumberOfPoints(); k++ ) {
393 MMG5_pPoint ppt = &mesh->point[k+1];
394 (*dataset)->GetPoint(k,mesh->point[k+1].c);
395 ppt->tag = MG_NUL;
396 ppt->tmp = 0;
397 ppt->ref = ptar->GetTuple1(k);
398 if ( ppt->ref < 0 ) {
399 ppt->ref = -ppt->ref;
400 ++nref;
401 }
402 }
403 }
404 else {
405 // read vertices only
406 for ( vtkIdType k = 0; k < (*dataset)->GetNumberOfPoints(); k++ ) {
407 MMG5_pPoint ppt = &mesh->point[k+1];
408 (*dataset)->GetPoint(k,mesh->point[k+1].c);
409 ppt->tag = MG_NUL;
410 ppt->tmp = 0;
411 ppt->ref = 0;
412 }
413 }
414
415 // Cells transfer toward Mmg data structure
416 assert ( !mesh->nei );
417 assert ( !mesh->nti );
418 assert ( !mesh->nai );
419 assert ( (mesh->npi == mesh->np) || !mesh->npi );
420
421 mesh->npi = 0;
422 MMG5_int nqi = 0;
423 MMG5_int npri = 0;
424 MMG5_int na = 0;
425 MMG5_int nbl_a = 0;
426 MMG5_int nt = 0;
427 MMG5_int nbl_t = 0;
428 MMG5_int ne = 0;
429
430 // Get pointer toward cells data containing element refs
431 vtkIdType numCells = (*dataset)->GetNumberOfCells();
432
433 if ( eltMeditRef > -1 ) {
434 // If present, read elements references
435 car = (*dataset)->GetCellData()->GetArray(eltMeditRef);
436
437 // Check the field name
438 assert ( strstr( car->GetName(), "medit:ref" ) );
439 // Check that we get 1 data only
440 assert ( car->GetNumberOfComponents() == 1 );
441
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);
446 return -1;
447 }
448 }
449
450 // Transfer cells to Mmg
451 for ( vtkIdType k = 0; k < numCells; k++ ) {
452 MMG5_pPoint ppt = NULL;
453 MMG5_pEdge pa = NULL;
454 MMG5_pTria ptt = NULL;
455 MMG5_pQuad pq = NULL;
456 MMG5_pTetra pt = NULL;
457 MMG5_pPrism ppr = NULL;
458
459 int typ = (*dataset)->GetCellType(k);
460 MMG5_int ref = 0;
461
462 switch ( typ ) {
463 case ( VTK_VERTEX ):
464 ppt = &mesh->point[++mesh->npi];
465 ppt->ref = car ? car->GetTuple1(k) : 0;
466 if ( ppt->ref < 0 ) {
467 ppt->ref = -ppt->ref;
468 ++nref;
469 }
470 break;
471
472 case ( VTK_POLY_LINE ):
473 int n;
474 n = (*dataset)->GetCell(k)->GetNumberOfPoints() - 1;
475
476 for ( int i=0; i<n; ++i ) {
477 ++mesh->nai;
478 ref = car ? car->GetTuple1(k) : 0;
479 }
480 /* Skip edges with iso ref */
481 if ( mesh->info.iso && MMG5_abs(ref) == mesh->info.isoref ) {
482 /* Skip this edge */
483 ++nbl_a;
484 }
485 else {
486 for ( int i=0; i < n; ++i ) {
487 pa = &mesh->edge[++na];
488 pa->a = (*dataset)->GetCell(k)->GetPointId(i) +1;
489 pa->b = (*dataset)->GetCell(k)->GetPointId(i+1)+1;
490 }
491 pa->ref = ref;
492 pa->tag |= MG_REF;
493 if ( pa->ref < 0 ) {
494 pa->ref = -pa->ref;
495 ++nref;
496 }
497 }
498 break;
499
500 case ( VTK_LINE ):
501 ++mesh->nai;
502 ref = car ? car->GetTuple1(k) : 0;
503
504 // Skip edges with iso ref
505 if ( mesh->info.iso && MMG5_abs(ref) == mesh->info.isoref ) {
506 /* Skip this edge */
507 ++nbl_a;
508 }
509 else {
510 pa = &mesh->edge[++na];
511 pa->a = (*dataset)->GetCell(k)->GetPointId(0)+1;
512 pa->b = (*dataset)->GetCell(k)->GetPointId(1)+1;
513 pa->ref = ref;
514 pa->tag |= MG_REF;
515 if ( pa->ref < 0 ) {
516 pa->ref = -pa->ref;
517 ++nref;
518 }
519 }
520 assert( na+nbl_a <= mesh->na );
521 break;
522
523 case ( VTK_TRIANGLE ):
524 ++mesh->nti;
525 ref = car ? car->GetTuple1(k) : 0;
526
527 // skip tria with iso ref in 3D
528 if ( mesh->info.iso && MMG5_abs(ref) == mesh->info.isoref && mesh->dim == 3 ) {
529 /* Skip this edge */
530 ++nbl_t;
531 }
532 else {
533 ptt = &mesh->tria[++nt];
534 for ( int i=0; i<3; ++i ) {
535 ptt->v[i] = (*dataset)->GetCell(k)->GetPointId(i)+1;
536 }
537 ptt->ref = ref;
538
539 if ( ptt->ref < 0 ) {
540 ptt->ref = -ptt->ref;
541 ++nref;
542 }
543 }
544 break;
545
546 case ( VTK_QUAD ):
547 pq = &mesh->quadra[++nqi];
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;
551 }
552 if ( pq->ref < 0 ) {
553 pq->ref = -pq->ref;
554 ++nref;
555 }
556 break;
557
558 case ( VTK_TETRA ):
559 // Skip volume elts if called from mmgs
560 if ( !mesh->tetra ) continue;
561
562 pt = &mesh->tetra[++mesh->nei];
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;
566 }
567 if ( pt->ref < 0 ) {
568 pt->ref = -pt->ref;
569 ++nref;
570 }
571 break;
572
573 case ( VTK_WEDGE ):
574 // Skip volume elts if called from mmgs
575 if ( !mesh->prism ) continue;
576
577 ppr = &mesh->prism[++npri];
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;
581 }
582 if ( ppr->ref < 0 ) {
583 ppr->ref = -ppr->ref;
584 ++nref;
585 }
586 break;
587
588 default:
589 if ( !mmgWarn1 ) {
590 printf(" ## Warning:%s: unexpected element type (%d).\n",
591 __func__,typ);
592 mmgWarn1 = 1;
593 }
594 }
595 }
596
597 // Checks
598 assert ( mesh->nai == mesh->na );
599 assert ( mesh->nti == mesh->nt );
600 assert ( mesh->nei == mesh->ne );
601 assert ( nqi == mesh->nquad );
602 assert ( npri == mesh->nprism );
603 assert ( na + nbl_a == mesh->na );
604 assert ( nt + nbl_t == mesh->nt );
605
606 // Array reallocation if ISO refs has been skipped
607 if ( mesh->info.iso ) {
608 if ( mesh->nt ) {
609 if( !nt )
611
612 else if ( nt < mesh->nt ) {
613 MMG5_ADD_MEM(mesh,(nt-mesh->nt)*sizeof(MMG5_Tria),"triangles",
614 fprintf(stderr," Exit program.\n");
615 return -1);
616 MMG5_SAFE_RECALLOC(mesh->tria,mesh->nt+1,(nt+1),MMG5_Tria,"triangles",
617 return -1);
618 }
619 mesh->nt = nt;
620 }
621 if ( mesh->na ) {
622 if( !na )
624 else if ( na < mesh->na ) {
625 MMG5_ADD_MEM(mesh,(na-mesh->na)*sizeof(MMG5_Edge),"edges",
626 fprintf(stderr," Exit program.\n");
627 return -1);
628 MMG5_SAFE_RECALLOC(mesh->edge,mesh->na+1,(na+1),MMG5_Edge,"edges",
629 return -1);
630 }
631 mesh->na = na;
632 }
633 }
634
636 if ( ier < 1 ) return ier;
637
638 if ( sol && sol[1] ) {
639 // Read the solution at nodes
640 // Init (*sol)[0] for the case where nsols=0
641 MMG5_pSol pmt = &sol[0][0];
642 MMG5_pSol psl = &sol[1][0];
643 psl->ver = mesh->ver;
644 psl->dim = mesh->dim;
645 psl->type = 1;
646
647 int isol = 0;
648 if ( nsols ) {
649 auto *pd = (*dataset)->GetPointData();
650 auto *cd = (*dataset)->GetCellData();
651
652 if ( pd ) {
653 int npointData = pd->GetNumberOfArrays();
654
655 for (int j = 0; j < npointData; j++) {
656 bool metricField = 0;
657 char chaine[MMG5_FILESTR_LGTH];
658 strcpy(chaine,pd->GetArrayName(j));
659
660 if ( strstr(chaine,"medit:ref" ) ) {
661 continue;
662 }
663 else if ( strstr(chaine,":metric") ) {
664 metricField = 1;
665 }
666
667 if ( metricField && (metricData*lsData) ) {
668 psl = sol[0];
669 }
670 else {
671 psl = &sol[1][isol];
672 ++isol;
673 }
674
675 psl->ver = mesh->ver;
676 psl->dim = mesh->dim;
677 psl->type = 1;
678 psl->entities = MMG5_Vertex;
679
680 if ( !MMG5_Set_inputSolName(mesh,psl,chaine) ) {
681 if ( !mmgWarn1 ) {
682 mmgWarn1 = 1;
683 fprintf(stderr,"\n ## Warning: %s: unable to set solution name for"
684 " at least 1 solution.\n",__func__);
685 }
686 }
687
688 if (pmt) {
689 if ( !MMG5_Set_inputSolName(mesh,pmt,chaine) ) {
690 if ( !mmgWarn1 ) {
691 mmgWarn1 = 1;
692 fprintf(stderr,"\n ## Warning: %s: unable to set solution name for"
693 " at least 1 solution.\n",__func__);
694 }
695 }
696 }
697
698 auto ar = pd->GetArray(j);
699
700 psl->np = ar->GetNumberOfTuples();
701 if ( mesh->np != psl->np ) {
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);
705 return -1;
706 }
707
708 int ncp = ar->GetNumberOfComponents();
709
710 if ( ncp == 1 ) {
711 psl->size = 1;
712 psl->type = 1;
713 }
714 else if ( ncp == 2 || ncp == 3 ) {
715 assert ( ncp == mesh->dim );
716 psl->size = ncp;
717 psl->type = 2;
718 }
719 else if ( ncp == (mesh->dim * mesh->dim) ) {
720 if ( metricData ) {
721 psl->size = (psl->dim*(psl->dim+1))/2;
722 psl->type = 3;
723 }
724 else {
725 psl->size = ncp;
726 psl->type = 4;
727 }
728 }
729 else {
730 fprintf(stderr," ** UNEXPECTED NUMBER OF COMPONENTS (%d). IGNORED \n",ncp);
731 return -1;
732 }
733
734 // mem alloc
735 if ( psl->m ) MMG5_DEL_MEM(mesh,psl->m);
736 psl->npmax = mesh->npmax;
737
738 MMG5_ADD_MEM(mesh,(psl->size*(psl->npmax+1))*sizeof(double),"initial solution",
739 fprintf(stderr," Exit program.\n");
740 return -1);
741 MMG5_SAFE_CALLOC(psl->m,psl->size*(psl->npmax+1),double,return -1);
742
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] );
748 }
749 break;
750
751 case ( 3 ):
752 // anisotropic sol
753 double dbuf[9];
754
755 for (MMG5_int k=1; k<=psl->np; k++) {
756 ar->GetTuple ( k-1, dbuf );
757 MMG5_int iadr = psl->size*k;
758
759 if ( !metricData ) {
760 // Non symmetric tensor
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];
766 }
767 else {
768 for ( int i=0 ; i<9 ; i++ ) {
769 psl->m[iadr+i] = dbuf[i];
770 }
771 }
772 }
773 else {
774 // Symmetric tensor
775 if ( psl->dim ==2 ) {
776 assert ( dbuf[1] == dbuf[2] );
777
778 psl->m[iadr] = dbuf[0];
779 psl->m[iadr+1] = dbuf[1];
780 psl->m[iadr+2] = dbuf[3];
781 }
782 else {
783 assert ( dbuf[1]==dbuf[3] || dbuf[2]==dbuf[6] || dbuf[5]==dbuf[7] );
784
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];
791 }
792 }
793 }
794
795 break;
796 default:
797 fprintf(stderr," ** UNEXPECTED METRIC TYPE (%d). EXIT PROGRAM \n",psl->type);
798 return -1;
799 }
800 }
801 }
802
803 if ( cd ) {
804 int ncellData = cd->GetNumberOfArrays();
805
806 for (int j = 0; j < ncellData; j++) {
807 // char *ptr = NULL;
808 char chaine[MMG5_FILESTR_LGTH];
809 strcpy(chaine,cd->GetArrayName(j));
810
811 if ( strstr(chaine,"medit:ref" ) ) {
812 continue;
813 }
814
815 psl = &sol[1][isol]; //*sol + isol;
816 psl->ver = mesh->ver;
817 psl->dim = mesh->dim;
818 psl->type = 1;
820
821 if ( !MMG5_Set_inputSolName(mesh,psl,chaine) ) {
822 if ( !mmgWarn1 ) {
823 mmgWarn1 = 1;
824 fprintf(stderr,"\n ## Warning: %s: unable to set solution name for"
825 " at least 1 solution.\n",__func__);
826 }
827 }
828
829 auto ar = cd->GetArray(j);
830
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);
836 return -1;
837 }
838
839 int ncp = ar->GetNumberOfComponents();
840
841 if ( ncp == 1 ) {
842 psl->size = 1;
843 psl->type = 1;
844 }
845 else if ( ncp == 2 || ncp == 3 ) {
846 assert ( ncp == mesh->dim );
847 psl->size = ncp;
848 psl->type = 2;
849 }
850 else if ( ncp == (mesh->dim * mesh->dim) ) {
851 psl->size = ncp;
852 psl->type = 3;
853 }
854 else {
855 fprintf(stderr," ** UNEXPECTED NUMBER OF COMPONENTS (%d). IGNORED \n",ncp);
856 return -1;
857 }
858
859 // mem alloc
860 if ( psl->m ) MMG5_DEL_MEM(mesh,psl->m);
861 psl->npmax = mesh->nemax;
862
863 MMG5_ADD_MEM(mesh,(psl->size*(psl->npmax+1))*sizeof(double),"initial solution",
864 fprintf(stderr," Exit program.\n");
865 return -1);
866 MMG5_SAFE_CALLOC(psl->m,psl->size*(psl->npmax+1),double,return -1);
867
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] );
873 }
874 break;
875
876 case ( 3 ):
877 // anisotropic sol
878 double dbuf[9];
879
880 for (MMG5_int k=1; k<=psl->np; k++) {
881 ar->GetTuple ( k-1, dbuf );
882 MMG5_int iadr = psl->size*k;
883
884 // Non symmetric tensor
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];
890 }
891 else {
892 for ( int i=0 ; i<9 ; i++ ) {
893 psl->m[iadr+i] = dbuf[i];
894 }
895 }
896 }
897
898 break;
899 default:
900 fprintf(stderr," ** UNEXPECTED METRIC TYPE (%d). EXIT PROGRAM \n",psl->type);
901 return -1;
902 }
903
904 ++isol;
905 }
906 }
907 }
908 }
909
910
911 (*dataset)->Delete();
912
913 return 1;
914}
int MMG5_Set_inputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solin)
int ier
MMG5_pMesh MMG5_pSol * sol
MMG5_pMesh * mesh
MMG5_pMesh char * filename
int MMG5_check_readedMesh(MMG5_pMesh mesh, MMG5_int nref)
Definition: inout.c:526
Types used throughout the Mmg libraries.
@ MMG5_Vertex
Definition: libmmgtypes.h:230
@ MMG5_Tetrahedron
Definition: libmmgtypes.h:233
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
#define MG_NUL
#define MMG5_ADD_MEM(mesh, size, message, law)
#define MMG5_FILESTR_LGTH
#define MG_REF
#define MMG5_DEL_MEM(mesh, ptr)
#define MMG5_SAFE_RECALLOC(ptr, prevSize, newSize, type, message, law)
Structure to store edges of am MMG mesh.
Definition: libmmgtypes.h:311
MMG5_int b
Definition: libmmgtypes.h:312
uint16_t tag
Definition: libmmgtypes.h:316
MMG5_int ref
Definition: libmmgtypes.h:313
MMG5_int a
Definition: libmmgtypes.h:312
int8_t iso
Definition: libmmgtypes.h:541
MMG5_int isoref
Definition: libmmgtypes.h:528
MMG mesh structure.
Definition: libmmgtypes.h:613
MMG5_pQuad quadra
Definition: libmmgtypes.h:656
MMG5_Info info
Definition: libmmgtypes.h:659
MMG5_int ne
Definition: libmmgtypes.h:620
MMG5_pPoint point
Definition: libmmgtypes.h:649
MMG5_pPrism prism
Definition: libmmgtypes.h:653
MMG5_int nquad
Definition: libmmgtypes.h:621
MMG5_int nemax
Definition: libmmgtypes.h:620
MMG5_int npmax
Definition: libmmgtypes.h:620
MMG5_int nei
Definition: libmmgtypes.h:620
MMG5_pTetra tetra
Definition: libmmgtypes.h:651
MMG5_int nt
Definition: libmmgtypes.h:620
MMG5_pTria tria
Definition: libmmgtypes.h:655
MMG5_int np
Definition: libmmgtypes.h:620
MMG5_pEdge edge
Definition: libmmgtypes.h:657
MMG5_int nti
Definition: libmmgtypes.h:620
MMG5_int npi
Definition: libmmgtypes.h:620
MMG5_int nai
Definition: libmmgtypes.h:620
MMG5_int nprism
Definition: libmmgtypes.h:621
MMG5_int na
Definition: libmmgtypes.h:620
Structure to store vertices of an MMG mesh.
Definition: libmmgtypes.h:276
MMG5_int tmp
Definition: libmmgtypes.h:286
double c[3]
Definition: libmmgtypes.h:277
uint16_t tag
Definition: libmmgtypes.h:290
MMG5_int ref
Definition: libmmgtypes.h:284
Structure to store prsim of a MMG mesh.
Definition: libmmgtypes.h:469
MMG5_int v[6]
Definition: libmmgtypes.h:470
MMG5_int ref
Definition: libmmgtypes.h:471
Structure to store quadrangles of an MMG mesh.
Definition: libmmgtypes.h:372
MMG5_int ref
Definition: libmmgtypes.h:374
MMG5_int v[4]
Definition: libmmgtypes.h:373
int entities
Definition: libmmgtypes.h:679
MMG5_int npmax
Definition: libmmgtypes.h:675
double * m
Definition: libmmgtypes.h:680
MMG5_int np
Definition: libmmgtypes.h:674
Structure to store tetrahedra of an MMG mesh.
Definition: libmmgtypes.h:407
MMG5_int v[4]
Definition: libmmgtypes.h:409
MMG5_int ref
Definition: libmmgtypes.h:410
Structure to store triangles of a MMG mesh.
Definition: libmmgtypes.h:338
MMG5_int ref
Definition: libmmgtypes.h:341
MMG5_int v[3]
Definition: libmmgtypes.h:340
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)
Definition: vtkparser.cpp:285
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)
Definition: vtkparser.cpp:245
vtkDataSet * MMG5_load_vtkXMLFile(const char *fileName)
I/O at VTK format.
Definition: vtkparser.cpp:63
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)
Definition: vtkparser.cpp:325
static int MMG5_count_vtkEntities(vtkDataSet *dataset, MMG5_pMesh mesh, int8_t *ptMeditRef, int8_t *eltMeditRef, int *nsols, int8_t *metricData, int8_t *lsData)
Definition: vtkparser.cpp:103
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)
Definition: vtkparser.cpp:366