Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
libmmg3d.c
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
42#include "libmmg3d.h"
45#include "mmgexterns_private.h"
46
51#define MMG5_RETURN_AND_PACK(mesh,met,sol,val)do \
52 { \
53 if ( !MMG3D_packMesh(mesh,met,sol) ) { \
54 mesh->npi = mesh->np; \
55 mesh->nti = mesh->nt; \
56 mesh->nai = mesh->na; \
57 mesh->nei = mesh->ne; \
58 mesh->xt = 0; \
59 if ( met ) { met->npi = met->np; } \
60 if ( sol ) { sol->npi = sol->np; } \
61 return MMG5_STRONGFAILURE; \
62 } \
63 _LIBMMG5_RETURN(mesh,met,sol,val); \
64 }while(0)
65
68 MMG5_int k;
69
70 mesh->xp = 0;
71 if ( mesh->adja )
73
75
76 if ( mesh->adjapr )
78
80
81 if ( mesh->xpoint )
83
84 for(k=1; k <=mesh->np; k++) {
85 mesh->point[k].xp = 0;
86 }
87
88 return;
89}
90
100 MMG5_pTetra pt;
101 MMG5_hgeom *ph;
102 MMG5_int k,nr;
103 int i;
104
105
106 /* rebuild triangles*/
107 if ( mesh->tria )
109 mesh->nt = 0;
110
111 if ( !MMG5_chkBdryTria(mesh) ) {
112 fprintf(stderr,"\n ## Error: %s: unable to rebuild triangles\n",__func__);
113 return -1;
114 }
115
116 /* build hash table for edges */
117 if ( mesh->htab.geom )
119
120 if ( mesh->edge )
122 mesh->na = 0;
123
124 nr = 0;
125 /* in the worst case (all edges are marked), we will have around 1 edge per *
126 * triangle (we count edges only one time) */
127 if ( MMG5_hNew(mesh,&mesh->htab,mesh->nt,3*(mesh->nt)) ) {
128 for (k=1; k<=mesh->ne; k++) {
129 pt = &mesh->tetra[k];
130 if ( MG_EOK(pt) && pt->xt ) {
131 for (i=0; i<6; i++) {
132 if ( mesh->xtetra[pt->xt].edg[i] ||
133 ( mesh->xtetra[pt->xt].tag[i] & MG_REQ ||
134 MG_EDG(mesh->xtetra[pt->xt].tag[i])) )
135 if ( !MMG5_hEdge(mesh,&mesh->htab,pt->v[MMG5_iare[i][0]],pt->v[MMG5_iare[i][1]],
136 mesh->xtetra[pt->xt].edg[i],mesh->xtetra[pt->xt].tag[i]))
137 return -1;
138 }
139 }
140 }
141
142 /* edges + ridges + required edges */
143 for (k=0; k<=mesh->htab.max; k++) {
144 ph = &mesh->htab.geom[k];
145 if ( !(ph->a) ) continue;
146 mesh->na++;
147 }
148 if ( mesh->na ) {
149 MMG5_ADD_MEM(mesh,(mesh->na+1)*sizeof(MMG5_Edge),"edges",
150 mesh->na = 0;
151 printf(" ## Warning: uncomplete mesh\n"));
152 }
153
154 if ( mesh->na ) {
155 MMG5_SAFE_CALLOC(mesh->edge,mesh->na+1,MMG5_Edge,return -1);
156
157 mesh->na = 0;
158 for (k=0; k<=mesh->htab.max; k++) {
159 ph = &mesh->htab.geom[k];
160 if ( !ph->a ) continue;
161 mesh->na++;
162 mesh->edge[mesh->na ].a = ph->a;
163 mesh->edge[mesh->na ].b = ph->b;
164 mesh->edge[mesh->na].tag = ( ph->tag | MG_REF ) ;
165 mesh->edge[mesh->na].ref = ph->ref;
166
167 if ( MG_GEO & ph->tag ) nr++;
168 }
169 }
171 }
172 else
173 mesh->memCur -= (long long)((3*mesh->nt+2)*sizeof(MMG5_hgeom));
174
175 mesh->nti = mesh->nt;
176 mesh->nai = mesh->na;
177
178 if ( mesh->info.imprim > 0 ) {
179 if ( mesh->na )
180 fprintf(stdout," NUMBER OF EDGES %8" MMG5_PRId " RIDGES %8" MMG5_PRId "\n",mesh->na,nr);
181 if ( mesh->nt )
182 fprintf(stdout," NUMBER OF TRIANGLES %8" MMG5_PRId "\n",mesh->nt);
183 }
184
185 return nr;
186}
187
188#ifndef USE_SCOTCH
200int MMG3D_mark_packedPoints(MMG5_pMesh mesh,MMG5_int *np,MMG5_int *nc) {
201 MMG5_pPoint ppt;
202 MMG5_int k;
203
204 for ( k=1; k<=mesh->np; ++k ) {
205 mesh->point[k].tmp = 0;
206 }
207
208 (*nc) = 0;
209 k = 1;
210 (*np) = mesh->np+1;
211
212 do {
213 ppt = &mesh->point[k];
214 if ( !MG_VOK(ppt) ) {
215 --(*np);
216 ppt = &mesh->point[*np];
217 assert ( ppt );
218
219 /* Search the last used point */
220 while ( !MG_VOK(ppt) && k < *np ) {
221 --(*np);
222 ppt = &mesh->point[*np];
223 }
224 }
225
226#ifndef NDEBUG
227 if ( !MG_VOK(ppt) ) {
228 assert ( k==*np );
229 }
230 if ( k==*np ) {
231 assert ( !MG_VOK(ppt) );
232 }
233#endif
234 if ( k==*np ) {
235 break;
236 }
237
238 ppt->tmp = k;
239
240 if ( ppt->tag & MG_NOSURF ) {
241 ppt->tag &= ~MG_NOSURF;
242 ppt->tag &= ~MG_REQ;
243 }
244
245 if ( ppt->tag & MG_CRN ) (*nc)++;
246
247 ppt->ref = MMG5_abs(ppt->ref);
248 }
249 while ( ++k < (*np) );
250
251 assert ( k==*np );
252 --(*np);
253
254#ifndef NDEBUG
255 for ( k=1; k<=mesh->np; ++k ) {
256 if ( MG_VOK(&mesh->point[k]) ) {
257 /* Check that all the used points have a tmp field */
258 assert(mesh->point[k].tmp);
259 }
260 }
261#endif
262
263 return 1;
264}
265
276 MMG5_pTetra pt,pt1;
277 MMG5_int k,iadr,iadr1,iadrv,*adjav,*adja,*adja1;
278 int i,voy;
279
280 if ( !mesh->ne ) {
281 return 1;
282 }
283
284 k = 1;
285 do {
286 pt = &mesh->tetra[k];
287 if ( !MG_EOK(pt) ) {
288 pt1 = &mesh->tetra[mesh->ne];
289 assert( pt && pt1 && MG_EOK(pt1) );
290 memcpy(pt,pt1,sizeof(MMG5_Tetra));
291
292 /* treat adja array */
293 iadr = 4*(k-1) + 1;
294 adja = &mesh->adja[iadr];
295 iadr1 = 4*(mesh->ne-1) + 1;
296 adja1 = &mesh->adja[iadr1];
297 for(i=0 ; i<4 ; i++) {
298 adja[i] = adja1[i];
299 if(!adja1[i]) continue;
300 iadrv = 4*(adja1[i]/4-1) + 1;
301 adjav = &mesh->adja[iadrv];
302 voy = i;
303 adjav[adja1[i]%4] = 4*k + voy;
304 }
305
306 if ( !MMG3D_delElt(mesh,mesh->ne) ) return 0;
307 }
308 }
309 while ( ++k < mesh->ne );
310
311 /* Recreate nil chain */
312 if ( mesh->ne >= mesh->nemax-1 )
313 mesh->nenil = 0;
314 else
315 mesh->nenil = mesh->ne + 1;
316
317 mesh->nei = mesh->ne;
318
319 if ( mesh->nenil )
320 for(k=mesh->nenil; k<mesh->nemax-1; k++)
321 mesh->tetra[k].v[3] = k+1;
322
323 return 1;
324}
325
336 MMG5_pTetra pt,pt1;
337 MMG5_int k;
338
339 if ( !mesh->ne ) {
340 return 1;
341 }
342
343 if ( mesh->tetra ) {
344 k = 1;
345 do {
346 pt = &mesh->tetra[k];
347 if ( !MG_EOK(pt) ) {
348 pt1 = &mesh->tetra[mesh->ne];
349 assert( pt && pt1 && MG_EOK(pt1) );
350 memcpy(pt,pt1,sizeof(MMG5_Tetra));
351 if ( !MMG3D_delElt(mesh,mesh->ne) ) return 0;
352 }
353 }
354 while ( ++k < mesh->ne );
355
356 /* Recreate nil chain */
357 if ( mesh->ne >= mesh->nemax-1 )
358 mesh->nenil = 0;
359 else
360 mesh->nenil = mesh->ne + 1;
361
362 if ( mesh->nenil )
363 for(k=mesh->nenil; k<mesh->nemax-1; k++)
364 mesh->tetra[k].v[0] = 0;
365 }
366 mesh->nei = mesh->ne;
367
368 return 1;
369}
370
380 MMG5_pPrism pp,pp1;
381 MMG5_pQuad pq,pq1;
382 MMG5_int k;
383
384 if ( mesh->prism ) {
385 k = 1;
386 do {
387 pp = &mesh->prism[k];
388 if ( !MG_EOK(pp) ) {
389 pp1 = &mesh->prism[mesh->nprism];
390 assert( pp && pp1 && MG_EOK(pp1) );
391 memcpy(pp,pp1,sizeof(MMG5_Prism));
392 --mesh->nprism;
393 }
394 }
395 while ( ++k < mesh->nprism );
396 }
397
398 if ( mesh->quadra ) {
399 k = 1;
400 do {
401 pq = &mesh->quadra[k];
402 if ( !MG_EOK(pq) ) {
403 pq1 = &mesh->quadra[mesh->nquad];
404 assert( pq && pq1 && MG_EOK(pq1) );
405 memcpy(pq,pq1,sizeof(MMG5_Quad));
406 --mesh->nquad;
407 }
408 }
409 while ( ++k < mesh->nquad );
410 }
411
412 return 1;
413}
414
424 MMG5_pPoint ppt,ppt1;
425 MMG5_int np,k,isol,isol1;
426 int i;
427
428 if ( sol && sol->m ) {
429 k = 1;
430 np = mesh->np;
431 /* Search the last used point */
432 while ( (!MG_VOK(mesh->point+np)) && np ) {
433 --np;
434 }
435
436 do {
437 ppt = &mesh->point[k];
438 if ( !MG_VOK(ppt) ) {
439 /* Copy last used point into first used one */
440 isol = k * sol->size;
441 isol1 = np * sol->size;
442
443 assert( sol->m+isol && sol->m+isol1 && MG_VOK(&mesh->point[np]) );
444
445 for (i=0; i<sol->size; i++) {
446 sol->m[isol + i] = sol->m[isol1 + i];
447 }
448
449 /* Search the last used point */
450 do {
451 --np;
452 ppt1 = &mesh->point[np];
453 }
454 while ( !MG_VOK(ppt1) && k < np );
455 }
456 }
457 while ( ++k < np );
458
459 sol->npi = sol->np = np;
460 }
461
462 return 1;
463}
464
473 MMG5_pPoint ppt,ppt1;
474 MMG5_int k;
475
476 k = 1;
477
478 while ( (!MG_VOK(&mesh->point[mesh->np])) && mesh->np ) {
480 }
481
482 do {
483 ppt = &mesh->point[k];
484 if ( !MG_VOK(ppt) ) {
485 ppt1 = &mesh->point[mesh->np];
486 assert( ppt && ppt1 && MG_VOK(ppt1) );
487 memcpy(ppt,ppt1,sizeof(MMG5_Point));
488
489 /* delete point without deleting xpoint */
490 memset(ppt1,0,sizeof(MMG5_Point));
491 ppt1->tag = MG_NUL;
492
493 while ( !MG_VOK((&mesh->point[mesh->np])) ) mesh->np--;
494
495 }
496
497 /* Copy the normal stored in the xpoint into ppt->n. */
498 if ( ppt->tag & MG_BDY &&
499 !(ppt->tag & MG_CRN || ppt->tag & MG_NOM || MG_EDG(ppt->tag)) ) {
500
501 if ( ppt->xp && mesh->xpoint ) {
502 memcpy(ppt->n,mesh->xpoint[ppt->xp].n1,3*sizeof(double));
503 ++mesh->nc1;
504 }
505 }
506 }
507 while ( ++k < mesh->np );
508
509 /* Recreate nil chain */
510 for(k=1 ; k<=mesh->np ; k++)
511 mesh->point[k].tmp = 0;
512
513 if ( mesh->np >= mesh->npmax-1 )
514 mesh->npnil = 0;
515 else
516 mesh->npnil = mesh->np + 1;
517
518 if ( mesh->npnil )
519 for(k=mesh->npnil; k<mesh->npmax-1; k++)
520 mesh->point[k].tmp = k+1;
521
522 /* Set mesh->npi to suitable value */
523 mesh->npi = mesh->np;
524
525 return 1;
526}
527
528#else
529
541int MMG3D_mark_packedPoints(MMG5_pMesh mesh,MMG5_int *np,MMG5_int *nc) {
542 MMG5_pPoint ppt;
543 MMG5_int k;
544
545 (*np) = (*nc) = 0;
546 for (k=1; k<=mesh->np; k++) {
547 ppt = &mesh->point[k];
548 if ( !MG_VOK(ppt) ) continue;
549 ppt->tmp = ++(*np);
550
551 if ( ppt->tag & MG_NOSURF ) {
552 ppt->tag &= ~MG_NOSURF;
553 ppt->tag &= ~MG_REQ;
554 }
555
556 if ( ppt->tag & MG_CRN ) (*nc)++;
557
558 ppt->ref = MMG5_abs(ppt->ref);
559 }
560 return 1;
561}
562
573 MMG5_pTetra pt,ptnew;
574 MMG5_int ne,nbl,k,iadr,iadrnew,iadrv,*adjav,*adja,*adjanew;
575 int i,voy;
576
577 ne = 0;
578 nbl = 1;
579 for (k=1; k<=mesh->ne; k++) {
580 pt = &mesh->tetra[k];
581 if ( !MG_EOK(pt) ) continue;
582
583 ne++;
584 if ( k!=nbl ) {
585 ptnew = &mesh->tetra[nbl];
586 memcpy(ptnew,pt,sizeof(MMG5_Tetra));
587
588 iadr = 4*(k-1) + 1;
589 adja = &mesh->adja[iadr];
590 iadrnew = 4*(nbl-1) + 1;
591 adjanew = &mesh->adja[iadrnew];
592 for(i=0 ; i<4 ; i++) {
593 adjanew[i] = adja[i];
594 if(!adja[i]) continue;
595 iadrv = 4*(adja[i]/4-1) +1;
596 adjav = &mesh->adja[iadrv];
597 voy = i;
598 adjav[adja[i]%4] = 4*nbl + voy;
599 }
600 }
601 nbl++;
602 }
603 mesh->nei = mesh->ne = ne;
604
605 /* Recreate nil chain */
606 if ( mesh->ne >= mesh->nemax-1 )
607 mesh->nenil = 0;
608 else
609 mesh->nenil = mesh->ne + 1;
610
611 if ( mesh->nenil )
612 for(k=mesh->nenil; k<mesh->nemax-1; k++)
613 mesh->tetra[k].v[3] = k+1;
614
615 return 1;
616}
617
628 MMG5_pTetra pt,ptnew;
629 MMG5_int ne,nbl,k;
630
631 ne = 0;
632 nbl = 1;
633 for (k=1; k<=mesh->ne; k++) {
634 pt = &mesh->tetra[k];
635 if ( !MG_EOK(pt) ) continue;
636
637 ne++;
638 if ( k!=nbl ) {
639 ptnew = &mesh->tetra[nbl];
640 memcpy(ptnew,pt,sizeof(MMG5_Tetra));
641 }
642 nbl++;
643 }
644 mesh->nei = mesh->ne = ne;
645
646 /* Recreate nil chain */
647 if ( mesh->ne >= mesh->nemax-1 )
648 mesh->nenil = 0;
649 else
650 mesh->nenil = mesh->ne + 1;
651
652 if ( mesh->nenil )
653 for(k=mesh->nenil; k<mesh->nemax-1; k++)
654 mesh->tetra[k].v[0] = 0;
655
656 return 1;
657}
658
668 MMG5_pPrism pp,ppnew;
669 MMG5_pQuad pq,pqnew;
670 MMG5_int k,ne,nbl;
671
672 ne = 0;
673 nbl = 1;
674 for (k=1; k<=mesh->nprism; k++) {
675 pp = &mesh->prism[k];
676 if ( !MG_EOK(pp) ) continue;
677
678 ++ne;
679 if ( k!=nbl ) {
680 ppnew = &mesh->prism[nbl];
681 memcpy(ppnew,pp,sizeof(MMG5_Prism));
682 }
683 ++nbl;
684 }
685 mesh->nprism = ne;
686
687 ne = 0;
688 nbl = 1;
689 for (k=1; k<=mesh->nquad; k++) {
690 pq = &mesh->quadra[k];
691 if ( !MG_EOK(pq) ) continue;
692
693 ++ne;
694 if ( k!=nbl ) {
695 pqnew = &mesh->quadra[nbl];
696 memcpy(pqnew,pq,sizeof(MMG5_Quad));
697 }
698 ++nbl;
699 }
700 mesh->nquad = ne;
701
702 return 1;
703}
704
714 MMG5_pPoint ppt;
715 MMG5_int k,isol,isolnew,np,nbl;
716 int i;
717
718 np = 0;
719 nbl = 1;
720 if ( sol && sol->m ) {
721 for (k=1; k<=mesh->np; k++) {
722 ppt = &mesh->point[k];
723 if ( !MG_VOK(ppt) ) continue;
724
725 ++np;
726
727 if ( k!= nbl ) {
728 isol = k * sol->size;
729 isolnew = nbl * sol->size;
730
731 for (i=0; i<sol->size; i++)
732 sol->m[isolnew + i] = sol->m[isol + i];
733 }
734 ++nbl;
735 }
736 sol->npi = sol->np = np;
737 }
738
739 return 1;
740}
741
750 MMG5_pPoint ppt,pptnew;
751 MMG5_int k,np,nbl;
752
753 nbl = 1;
754 mesh->nc1 = 0;
755 np = 0;
756 for (k=1; k<=mesh->np; k++) {
757 ppt = &mesh->point[k];
758 if ( !MG_VOK(ppt) ) continue;
759
760 if ( ppt->tag & MG_BDY &&
761 !(ppt->tag & MG_CRN || ppt->tag & MG_NOM || MG_EDG(ppt->tag)) ) {
762
763 if ( ppt->xp && mesh->xpoint ) {
764 memcpy(ppt->n,mesh->xpoint[ppt->xp].n1,3*sizeof(double));
765 ++mesh->nc1;
766 }
767 }
768
769 np++;
770 if ( k!=nbl ) {
771 pptnew = &mesh->point[nbl];
772 memmove(pptnew,ppt,sizeof(MMG5_Point));
773 memset(ppt,0,sizeof(MMG5_Point));
774 ppt->tag = MG_NUL;
775 }
776
777 nbl++;
778 }
779 mesh->npi = mesh->np = np;
780
781 for(k=1 ; k<=mesh->np ; k++)
782 mesh->point[k].tmp = 0;
783
784 if ( mesh->np >= mesh->npmax-1 )
785 mesh->npnil = 0;
786 else
787 mesh->npnil = mesh->np + 1;
788
789 if ( mesh->npnil )
790 for(k=mesh->npnil; k<mesh->npmax-1; k++)
791 mesh->point[k].tmp = k+1;
792
793 return 1;
794}
795
796#endif
797
798
808 MMG5_pTetra pt;
809 MMG5_pPrism pp;
810 MMG5_pQuad pq;
811 MMG5_int k;
812
813 for (k=1; k<=mesh->ne; k++) {
814 pt = &mesh->tetra[k];
815 if ( !MG_EOK(pt) ) continue;
816
817 pt->v[0] = mesh->point[pt->v[0]].tmp;
818 pt->v[1] = mesh->point[pt->v[1]].tmp;
819 pt->v[2] = mesh->point[pt->v[2]].tmp;
820 pt->v[3] = mesh->point[pt->v[3]].tmp;
821 }
822 for (k=1; k<=mesh->nprism; k++) {
823 pp = &mesh->prism[k];
824 if ( !MG_EOK(pp) ) continue;
825
826 pp->v[0] = mesh->point[pp->v[0]].tmp;
827 pp->v[1] = mesh->point[pp->v[1]].tmp;
828 pp->v[2] = mesh->point[pp->v[2]].tmp;
829 pp->v[3] = mesh->point[pp->v[3]].tmp;
830 pp->v[4] = mesh->point[pp->v[4]].tmp;
831 pp->v[5] = mesh->point[pp->v[5]].tmp;
832 }
833 for (k=1; k<=mesh->nquad; k++) {
834 pq = &mesh->quadra[k];
835 if ( !MG_EOK(pq) ) continue;
836
837 pq->v[0] = mesh->point[pq->v[0]].tmp;
838 pq->v[1] = mesh->point[pq->v[1]].tmp;
839 pq->v[2] = mesh->point[pq->v[2]].tmp;
840 pq->v[3] = mesh->point[pq->v[3]].tmp;
841 }
842 return 1;
843}
844
854 MMG5_int np, nc;
855
857 if ( !MMG3D_mark_packedPoints(mesh,&np,&nc) ) return -1;
858
860 if ( !MMG3D_update_eltsVertices(mesh) ) return -1;
861
863 if ( MMG3D_pack_pointArray(mesh)<0 ) return -1;
864
865 return nc;
866}
867
876 MMG5_pTetra pt;
877 MMG5_int k;
878 int i;
879
880 for (k=1; k<=mesh->ne; k++) {
881 pt = &mesh->tetra[k];
882 if ( MG_EOK(pt) && pt->xt ) {
883
884 for (i=0; i<6; i++) {
885 if ( mesh->xtetra[pt->xt].tag[i] & MG_NOSURF ) {
886 mesh->xtetra[pt->xt].tag[i] &= ~MG_REQ;
887 mesh->xtetra[pt->xt].tag[i] &= ~MG_NOSURF;
888 }
889 }
890 }
891 }
892 return;
893}
894
906 MMG5_int nc,nr;
907
908 /* Remove non wanted subdomains if needed */
910
911 /* compact vertices */
912 if ( !mesh->point ) {
913 fprintf(stderr, "\n ## Error: %s: points array not allocated.\n",
914 __func__);
915 return 0;
916 }
917 if ( !mesh->tetra ) {
918 fprintf(stderr, "\n ## Error: %s: tetra array not allocated.\n",
919 __func__);
920 return 0;
921 }
922
923 /* compact tetrahedra */
924 if ( mesh->adja ) {
925 if ( !MMG3D_pack_tetraAndAdja(mesh) ) return 0;
926 }
927 else {
928 if ( !MMG3D_pack_tetra(mesh) ) return 0;
929 }
930
931 /* update prisms and quads vertex indices */
932 if ( !MMG3D_pack_prismsAndQuads(mesh) ) return 0;
933
934 /* compact solutions (metric, level-set, displacement...) */
935 if ( met && met->m )
936 if ( !MMG3D_pack_sol(mesh,met) ) return 0;
937
938 if ( sol && sol->m )
939 if ( !MMG3D_pack_sol(mesh,sol) ) return 0;
940
941 /*compact vertices*/
943 if ( nc<0 ) return 0;
944
945 if ( met && met->m ) assert(met->np ==mesh->np);
946 if ( sol && sol->m ) assert(sol->np == mesh->np);
947
948 /* create prism adjacency */
949 if ( !MMG3D_hashPrism(mesh) ) {
950 fprintf(stderr,"\n ## Error: %s: prism hashing problem. Exit program.\n",
951 __func__);
952 return 0;
953 }
954
955 /* Remove the MG_REQ tags added by the nosurf option */
957
958 if ( mesh->info.imprim > 0 ) {
959 fprintf(stdout," NUMBER OF VERTICES %8" MMG5_PRId " CORNERS %8" MMG5_PRId "\n",mesh->np,nc);
960 fprintf(stdout," NUMBER OF TETRAHEDRA %8" MMG5_PRId "\n",mesh->ne);
961 }
962
963 nr = MMG3D_bdryBuild(mesh);
964 if ( nr < 0 ) return 0;
965
966 /* to could save the mesh, the adjacency have to be correct */
967 if ( mesh->info.ddebug && (!MMG5_chkmsh(mesh,1,1) ) ) {
968 fprintf(stderr,"\n ## Warning: %s: invalid mesh.\n",__func__);
969 return 0;
970 }
971
972 return 1;
973}
974
976 MMG5_pSol sol=NULL; // unused
977 mytime ctim[TIMEMAX];
978 char stim[32];
979
981 assert ( mesh );
982 assert ( met );
983 assert ( mesh->point );
984 assert ( mesh->tetra );
985
986 MMG5_version(mesh,"3D");
987
989
990
992
996
997 signal(SIGABRT,MMG5_excfun);
998 signal(SIGFPE,MMG5_excfun);
999 signal(SIGILL,MMG5_excfun);
1000 signal(SIGSEGV,MMG5_excfun);
1001 signal(SIGTERM,MMG5_excfun);
1002 signal(SIGINT,MMG5_excfun);
1003
1004 tminit(ctim,TIMEMAX);
1005 chrono(ON,&(ctim[0]));
1006
1007 /* Check options */
1008 if ( mesh->info.lag > -1 ) {
1009 fprintf(stderr,"\n ## ERROR: LAGRANGIAN MODE UNAVAILABLE (MMG3D_IPARAM_lag):\n"
1010 " YOU MUST CALL THE MMG3D_MMG3DMOV FUNCTION TO MOVE A RIGIDBODY.\n");
1012 }
1013 else if ( mesh->info.iso || mesh->info.isosurf ) {
1014 fprintf(stderr,"\n ## ERROR: LEVEL-SET DISCRETISATION UNAVAILABLE"
1015 " (MMG3D_IPARAM_iso or MMG3D_IARAM_isosurf ):\n"
1016 " YOU MUST CALL THE MMG3D_MMG3DMOV FUNCTION TO USE THIS OPTION.\n");
1018 }
1019 else if ( mesh->info.optimLES && met->size==6 ) {
1020 fprintf(stdout,"\n ## ERROR: STRONG MESH OPTIMIZATION FOR LES METHODS"
1021 " UNAVAILABLE (MMG3D_IPARAM_optimLES) WITH AN ANISOTROPIC METRIC.\n");
1023 }
1024
1025 if ( mesh->info.imprim > 0 ) fprintf(stdout,"\n -- MMG3DLIB: INPUT DATA\n");
1026
1027 chrono(ON,&(ctim[1]));
1028
1029 /* check input */
1030 if ( met->np && (met->np != mesh->np) ) {
1031 fprintf(stdout,"\n ## WARNING: WRONG SOLUTION NUMBER. IGNORED\n");
1032 MMG5_DEL_MEM(mesh,met->m);
1033 met->np = 0;
1034 }
1035 else if ( met->size!=1 && met->size!=6 ) {
1036 fprintf(stderr,"\n ## ERROR: WRONG DATA TYPE.\n");
1038 }
1039
1040 /* specific meshing */
1041 if ( met->np ) {
1042 if ( mesh->info.optim ) {
1043 printf("\n ## ERROR: MISMATCH OPTIONS: OPTIM OPTION CAN NOT BE USED"
1044 " WITH AN INPUT METRIC.\n");
1046 }
1047
1048 if ( mesh->info.hsiz>0. ) {
1049 printf("\n ## ERROR: MISMATCH OPTIONS: HSIZ OPTION CAN NOT BE USED"
1050 " WITH AN INPUT METRIC.\n");
1052 }
1053 }
1054
1055 if ( mesh->info.optim && mesh->info.hsiz>0. ) {
1056 printf("\n ## ERROR: MISMATCH OPTIONS: HSIZ AND OPTIM OPTIONS CAN NOT BE USED"
1057 " TOGETHER.\n");
1059 }
1060
1061#ifdef USE_SCOTCH
1063#endif
1064
1065 chrono(OFF,&(ctim[1]));
1066 printim(ctim[1].gdif,stim);
1067 if ( mesh->info.imprim > 0 )
1068 fprintf(stdout," -- INPUT DATA COMPLETED. %s\n",stim);
1069
1070 /* analysis */
1071 chrono(ON,&(ctim[2]));
1072 if ( mesh->info.imprim > 0 ) {
1073 fprintf(stdout,"\n -- PHASE 1 : ANALYSIS\n");
1074 }
1075
1076 /* reset fem value to user setting (needed for multiple library call) */
1078
1079 /* scaling mesh */
1081
1082 MMG3D_setfunc(mesh,met);
1083
1084 /* specific meshing */
1085 if ( mesh->info.optim ) {
1086 if ( !MMG3D_doSol(mesh,met) ) {
1089 }
1090 }
1091
1092 if ( mesh->info.hsiz > 0. ) {
1093 if ( !MMG3D_Set_constantSize(mesh,met) ) {
1096 }
1097 }
1098
1100
1101 if ( mesh->info.imprim > 0 || mesh->info.imprim < -1 ) {
1102 if ( !MMG3D_inqua(mesh,met) ) {
1105 }
1106 }
1107
1108 /* mesh analysis */
1109 if ( !MMG3D_analys(mesh) ) {
1112 }
1113
1114 if ( mesh->info.imprim > 1 && met->m ) MMG3D_prilen(mesh,met,0);
1115
1116 chrono(OFF,&(ctim[2]));
1117 printim(ctim[2].gdif,stim);
1118 if ( mesh->info.imprim > 0 )
1119 fprintf(stdout," -- PHASE 1 COMPLETED. %s\n",stim);
1120
1121 /* mesh adaptation */
1122 chrono(ON,&(ctim[3]));
1123 if ( mesh->info.imprim > 0 ) {
1124 fprintf(stdout,"\n -- PHASE 2 : %s MESHING\n",met->size < 6 ? "ISOTROPIC" : "ANISOTROPIC");
1125 }
1126
1127 /* renumerotation if available */
1128 if ( !MMG5_scotchCall(mesh,met,NULL,NULL) )
1129 {
1132 }
1133
1134#ifdef MMG_PATTERN
1135 if ( !MMG5_mmg3d1_pattern(mesh,met,NULL) ) {
1136 if ( !(mesh->adja) && !MMG3D_hashTetra(mesh,1) ) {
1137 fprintf(stderr,"\n ## Hashing problem. Invalid mesh.\n");
1139 }
1142 }
1143#else
1144 if ( !MMG5_mmg3d1_delone(mesh,met,NULL) ) {
1145 if ( (!mesh->adja) && !MMG3D_hashTetra(mesh,1) ) {
1146 fprintf(stderr,"\n ## Hashing problem. Invalid mesh.\n");
1148 }
1151 }
1152#endif
1153
1154 chrono(OFF,&(ctim[3]));
1155 printim(ctim[3].gdif,stim);
1156 if ( mesh->info.imprim > 0 ) {
1157 fprintf(stdout," -- PHASE 2 COMPLETED. %s\n",stim);
1158 }
1159
1160 /* last renum to give back a good numbering to the user */
1161 if ( !MMG5_scotchCall(mesh,met,NULL,NULL) )
1162 {
1165 }
1166
1167 /* save file */
1168 if ( !MMG3D_outqua(mesh,met) ) {
1171 }
1172
1173 if ( mesh->info.imprim > 4 && met->m )
1174 MMG3D_prilen(mesh,met,1);
1175
1176 chrono(ON,&(ctim[1]));
1177 if ( mesh->info.imprim > 0 ) fprintf(stdout,"\n -- MESH PACKED UP\n");
1180 chrono(OFF,&(ctim[1]));
1181
1182 chrono(OFF,&ctim[0]);
1183 printim(ctim[0].gdif,stim);
1184 if ( mesh->info.imprim >= 0 ) {
1185 fprintf(stdout,"\n MMG3DLIB: ELAPSED TIME %s\n",stim);
1186 fprintf(stdout,"\n %s\n END OF MODULE MMG3D\n %s\n\n",MG_STR,MG_STR);
1187 }
1189}
1190
1191/* Level set function discretization and remeshing mode */
1193 MMG5_pSol met=NULL;
1194 mytime ctim[TIMEMAX];
1195 char stim[32];
1196 int8_t mettofree = 0;
1197
1199 assert ( mesh );
1200 assert ( sol );
1201 assert ( mesh->point );
1202 assert ( mesh->tetra );
1203
1204 if ( (!mesh->info.iso) && (!mesh->info.isosurf) ) {
1205 fprintf(stdout,"\n ## WARNING: ISO MODE NOT PROVIDED: ENABLING ISOVALUE DISCRETIZATION MODE (-ls) \n");
1206 mesh->info.iso = 1;
1207 }
1208
1209 if ( !umet ) {
1210 /* User doesn't provide the metric (library mode only), allocate our own one */
1212 mettofree = 1;
1213 }
1214 else {
1215 /* Using appli we always pass here */
1216 met = umet;
1217 }
1218
1220
1221 signal(SIGABRT,MMG5_excfun);
1222 signal(SIGFPE,MMG5_excfun);
1223 signal(SIGILL,MMG5_excfun);
1224 signal(SIGSEGV,MMG5_excfun);
1225 signal(SIGTERM,MMG5_excfun);
1226 signal(SIGINT,MMG5_excfun);
1227
1228 tminit(ctim,TIMEMAX);
1229 chrono(ON,&(ctim[0]));
1230
1231 /* Check options */
1232 if ( mesh->info.lag > -1 ) {
1233 fprintf(stderr,"\n ## ERROR: LAGRANGIAN MODE UNAVAILABLE (MMG3D_IPARAM_lag):\n"
1234 " YOU MUST CALL THE MMG3D_MMG3DMOV FUNCTION TO MOVE A RIGIDBODY.\n");
1235 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1237 }
1238 else if ( mesh->info.optimLES ) {
1239 fprintf(stdout,"\n ## ERROR: STRONG MESH OPTIMIZATION FOR LES METHODS"
1240 " UNAVAILABLE (MMG3D_IPARAM_optimLES) IN ISOSURFACE"
1241 " DISCRETIZATION MODE.\n");
1242 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1244 }
1245
1246 /* Specific meshing */
1247 if ( met && met->np ) {
1248 if ( mesh->info.optim ) {
1249 printf("\n ## ERROR: MISMATCH OPTIONS: OPTIM OPTION CAN NOT BE USED"
1250 " WITH AN INPUT METRIC.\n");
1251 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1253 }
1254
1255 if ( mesh->info.hsiz>0. ) {
1256 printf("\n ## ERROR: MISMATCH OPTIONS: HSIZ OPTION CAN NOT BE USED"
1257 " WITH AN INPUT METRIC.\n");
1258 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1260 }
1261 }
1262
1263 if ( mesh->info.optim && mesh->info.hsiz>0. ) {
1264 printf("\n ## ERROR: MISMATCH OPTIONS: HSIZ AND OPTIM OPTIONS CAN NOT BE USED"
1265 " TOGETHER.\n");
1266 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1268 }
1269
1270#ifdef USE_SCOTCH
1272#endif
1273
1274 if ( mesh->info.imprim > 0 ) fprintf(stdout,"\n -- MMG3DLS: INPUT DATA\n");
1275 /* load data */
1276 chrono(ON,&(ctim[1]));
1278
1282
1283 if ( sol->np && (sol->np != mesh->np) ) {
1284 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1286 }
1287 else if ( sol->size!=1 ) {
1288 fprintf(stderr,"\n ## ERROR: WRONG DATA TYPE.\n");
1289 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1291 }
1292 if ( met && met->np && (met->np != mesh->np) ) {
1293 fprintf(stdout,"\n ## WARNING: WRONG METRIC NUMBER. IGNORED\n");
1294 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1296 }
1297
1298 /* clean old isosurface */
1299 if ( !MMG3D_Clean_isoSurf(mesh) ) {
1300 fprintf(stderr,"\n ## Unable to clean old isosurface.\n");
1302 }
1303
1304 chrono(OFF,&(ctim[1]));
1305 printim(ctim[1].gdif,stim);
1306 if ( mesh->info.imprim > 0 )
1307 fprintf(stdout," -- INPUT DATA COMPLETED. %s\n",stim);
1308
1309 chrono(ON,&(ctim[2]));
1310
1311 if ( mesh->info.imprim > 0 ) {
1312 fprintf(stdout,"\n -- PHASE 1 : ISOSURFACE DISCRETIZATION\n");
1313 }
1314
1315 /* reset fem value to user setting (needed for multiple library call) */
1317
1318 /* scaling mesh */
1319 if ( !MMG5_scaleMesh(mesh,met,sol) ) {
1320 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1322 }
1323
1324 MMG3D_setfunc(mesh,met);
1325
1326 if ( !MMG3D_tetraQual(mesh,met,0) ) {
1327 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1329 }
1330
1331 if ( mesh->info.imprim > 0 || mesh->info.imprim < -1 ) {
1332 if ( !MMG3D_inqua(mesh,met) ) {
1333 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1334 if ( !MMG5_unscaleMesh(mesh,met,sol) ) {
1336 }
1338 }
1339 }
1340
1341 if ( !sol->np ) {
1342 fprintf(stderr,"\n ## ERROR: A VALID SOLUTION FILE IS NEEDED \n");
1343 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1345 }
1346
1347 /* Specific meshing: compute optim option here because after isovalue
1348 * discretization mesh elements have too bad qualities */
1349 if ( mesh->info.optim ) {
1350 if ( !MMG3D_doSol(mesh,met) ) {
1351 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1352 if ( !MMG5_unscaleMesh(mesh,met,sol) ) {
1355 }
1356 }
1357
1358 /* Discretization of the mesh->info.ls isovalue of sol in the mesh */
1359 if ( !MMG3D_mmg3d2(mesh,sol,met) ) {
1360 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1362 }
1363
1364 chrono(OFF,&(ctim[2]));
1365 printim(ctim[2].gdif,stim);
1366 if ( mesh->info.imprim > 0 )
1367 fprintf(stdout," -- PHASE 1 COMPLETED. %s\n",stim);
1368
1369 chrono(ON,&(ctim[3]));
1370 if ( mesh->info.imprim > 0 ) {
1371 fprintf(stdout,"\n -- PHASE 2 : ANALYSIS\n");
1372 }
1373
1374 /* Specific meshing: compute hsiz on mesh after isovalue discretization */
1375 if ( mesh->info.hsiz > 0. ) {
1376 if ( !MMG3D_Set_constantSize(mesh,met) ) {
1377 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1380 }
1381 }
1382
1383 /* mesh analysis */
1384 if ( !MMG3D_analys(mesh) ) {
1385 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1388 }
1389
1390
1391 chrono(OFF,&(ctim[3]));
1392 printim(ctim[3].gdif,stim);
1393 if ( mesh->info.imprim > 0 )
1394 fprintf(stdout," -- PHASE 2 COMPLETED. %s\n",stim);
1395
1396 /* mesh adaptation */
1397 chrono(ON,&(ctim[4]));
1398 if ( mesh->info.imprim > 0 ) {
1399 fprintf(stdout,"\n -- PHASE 3 : MESH IMPROVEMENT\n");
1400 }
1401
1402 /* renumerotation if available */
1403 if ( !MMG5_scotchCall(mesh,met,NULL,NULL) )
1404 {
1405 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1408 }
1409
1410#ifdef MMG_PATTERN
1411 if ( !MMG5_mmg3d1_pattern(mesh,met,NULL) ) {
1412 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1413 if ( !(mesh->adja) && !MMG3D_hashTetra(mesh,1) ) {
1414 fprintf(stderr,"\n ## Hashing problem. Invalid mesh.\n");
1416 }
1419 }
1420#else
1421 if ( !MMG5_mmg3d1_pattern(mesh,met,NULL) ) {
1422 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1423 if ( !(mesh->adja) && !MMG3D_hashTetra(mesh,1) ) {
1424 fprintf(stderr,"\n ## Hashing problem. Invalid mesh.\n");
1426 }
1429 }
1430#endif
1431
1432 chrono(OFF,&(ctim[4]));
1433 printim(ctim[4].gdif,stim);
1434 if ( mesh->info.imprim > 0 ) {
1435 fprintf(stdout," -- PHASE 3 COMPLETED. %s\n",stim);
1436 }
1437
1438 /* last renum to give back a good numbering to the user */
1439 if ( !MMG5_scotchCall(mesh,met,NULL,NULL) )
1440 {
1443 }
1444
1445 /* save file */
1446 if ( !MMG3D_outqua(mesh,met) ) {
1447 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1450 }
1451
1452 chrono(ON,&(ctim[1]));
1453 if ( mesh->info.imprim > 0 ) fprintf(stdout,"\n -- MESH PACKED UP\n");
1454 if ( !MMG5_unscaleMesh(mesh,met,sol) ) {
1455 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1457 }
1458 if ( !MMG3D_packMesh(mesh,sol,met) ) {
1459 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1461 }
1462 chrono(OFF,&(ctim[1]));
1463
1464 chrono(OFF,&ctim[0]);
1465 printim(ctim[0].gdif,stim);
1466 if ( mesh->info.imprim >= 0 ) {
1467 fprintf(stdout,"\n MMG3DLS: ELAPSED TIME %s\n",stim);
1468 fprintf(stdout,"\n %s\n END OF MODULE MMG3D\n %s\n\n",MG_STR,MG_STR);
1469 }
1470 if ( mettofree ) { MMG5_DEL_MEM(mesh,met->m);MMG5_SAFE_FREE (met); }
1472}
1473
1474
1476 mytime ctim[TIMEMAX];
1477 char stim[32];
1478 MMG5_int *invalidTets;
1479 MMG5_int k;
1480 MMG5_int ier;
1481
1483 assert ( mesh );
1484 assert ( met );
1485 assert ( mesh->point );
1486 assert ( mesh->tetra );
1487
1488 MMG5_version(mesh,"3D");
1489
1491
1492 signal(SIGABRT,MMG5_excfun);
1493 signal(SIGFPE,MMG5_excfun);
1494 signal(SIGILL,MMG5_excfun);
1495 signal(SIGSEGV,MMG5_excfun);
1496 signal(SIGTERM,MMG5_excfun);
1497 signal(SIGINT,MMG5_excfun);
1498
1499 tminit(ctim,TIMEMAX);
1500 chrono(ON,&(ctim[0]));
1501
1502 /* Check options */
1503 if ( mesh->info.iso || mesh->info.isosurf ) {
1504 fprintf(stderr,"\n ## ERROR: LEVEL-SET DISCRETISATION UNAVAILABLE"
1505 " (MMG3D_IPARAM_iso || MMG3D_IPARAM_isosurf ):\n"
1506 " YOU MUST CALL THE MMG3D_mmg3dls FUNCTION TO USE THIS OPTION.\n");
1508 }
1509 else if ( mesh->info.optimLES ) {
1510 fprintf(stdout,"\n ## ERROR: STRONG MESH OPTIMIZATION FOR LES METHODS"
1511 " UNAVAILABLE (MMG3D_IPARAM_optimLES) IN LAGRANGIAN MODE.\n");
1513 }
1514 if ( mesh->info.optim ) {
1515 printf("\n ## ERROR: OPTIM OPTION UNAVAILABLE IN LAGRANGIAN MODE\n");
1517 }
1518 if ( mesh->info.hsiz>0. ) {
1519 printf("\n ## ERROR: HSIZ OPTION UNAVAILABLE IN LAGRANGIAN MODE\n");
1521 }
1522
1523#ifdef USE_SCOTCH
1525#endif
1526
1527 if ( mesh->info.imprim > 0 ) fprintf(stdout,"\n -- MMG3DMOV: INPUT DATA\n");
1528 /* load data */
1529 chrono(ON,&(ctim[1]));
1531
1535
1536 if ( mesh->info.lag == -1 ) {
1537 if ( mesh->info.imprim > 0 )
1538 fprintf(stdout,"\n ## Warning: displacement mode for the rigidbody"
1539 " movement is not set.\n"
1540 " Lagrangian displacement computed according"
1541 " to mode 1.\n");
1542 mesh->info.lag = 1;
1543 }
1544
1545#ifndef USE_ELAS
1546 fprintf(stderr,"\n ## ERROR: YOU NEED TO COMPILE WITH THE USE_ELAS"
1547 " CMake's FLAG SET TO ON TO USE THE RIGIDBODY MOVEMENT LIBRARY.\n");
1549#endif
1550
1551 if ( !disp ) {
1552 fprintf(stderr,"\n ## ERROR: IN LAGRANGIAN MODE, A STRUCTURE OF TYPE"
1553 " \"MMG5_pSoL\" IS NEEDED TO STORE THE DISPLACEMENT FIELD.\n"
1554 " THIS STRUCTURE MUST BE DIFFERENT FROM THE ONE USED"
1555 " TO STORE THE METRIC.\n");
1557 }
1558 if (disp->np && (disp->np != mesh->np) ) {
1559 fprintf(stdout,"\n ## WARNING: WRONG SOLUTION NUMBER. IGNORED\n");
1560 MMG5_DEL_MEM(mesh,disp->m);
1561 disp->np = 0;
1562 }
1563 else if (disp->size!=3) {
1564 fprintf(stderr,"\n ## ERROR: LAGRANGIAN MOTION OPTION NEED A VECTORIAL DISPLACEMENT\n");
1566 }
1567
1568 chrono(OFF,&(ctim[1]));
1569 printim(ctim[1].gdif,stim);
1570 if ( mesh->info.imprim > 0 )
1571 fprintf(stdout," -- INPUT DATA COMPLETED. %s\n",stim);
1572
1573 /* analysis */
1574 chrono(ON,&(ctim[2]));
1575
1576 if ( mesh->info.imprim > 0 ) {
1577 fprintf(stdout,"\n -- PHASE 1 : ANALYSIS\n");
1578 }
1579
1580 /* reset fem value to user setting (needed for multiple library call) */
1582
1583 /* scaling mesh */
1584 if ( !MMG5_scaleMesh(mesh,met,disp) ) _LIBMMG5_RETURN(mesh,met,disp,MMG5_STRONGFAILURE);
1585
1586 MMG3D_setfunc(mesh,met);
1587
1588 if ( !MMG3D_tetraQual(mesh,met,0) ) _LIBMMG5_RETURN(mesh,met,disp,MMG5_LOWFAILURE);
1589
1590 if ( mesh->info.imprim > 0 || mesh->info.imprim < -1 ) {
1591 if ( !MMG3D_inqua(mesh,met) ) {
1593 }
1594 }
1595
1596 /* mesh analysis */
1597 if ( !MMG3D_analys(mesh) ) {
1598 if ( !MMG5_unscaleMesh(mesh,met,disp) ) _LIBMMG5_RETURN(mesh,met,disp,MMG5_STRONGFAILURE);
1600 }
1601
1602 if ( mesh->info.imprim > 4 && met->m ) {
1603 MMG3D_prilen(mesh,met,0);
1604 }
1605
1606 chrono(OFF,&(ctim[2]));
1607 printim(ctim[2].gdif,stim);
1608 if ( mesh->info.imprim > 0 )
1609 fprintf(stdout," -- PHASE 1 COMPLETED. %s\n",stim);
1610
1611 /* mesh adaptation */
1612 chrono(ON,&(ctim[3]));
1613 if ( mesh->info.imprim > 0 ) {
1614 fprintf(stdout,"\n -- PHASE 2 : LAGRANGIAN MOTION\n");
1615 }
1616
1617 /* renumerotation if available */
1618 if ( !MMG5_scotchCall(mesh,met,disp,NULL) )
1619 {
1620 if ( !MMG5_unscaleMesh(mesh,met,disp) ) _LIBMMG5_RETURN(mesh,met,disp,MMG5_STRONGFAILURE);
1622 }
1623
1624 /* Lagrangian mode */
1625 invalidTets = NULL;
1626 ier = MMG5_mmg3d3(mesh,disp,met,&invalidTets);
1627 if ( !ier ) {
1628 disp->npi = disp->np;
1630 }
1631 else if ( ier < 0 ) {
1632 printf("\n ## Warning: Unable to perform any movement "
1633 "(%" MMG5_PRId " intersecting tetrahedra).\n",-ier);
1634 if ( mesh->info.imprim > 1 ) {
1635 printf(" List of invalid tets: ");
1636 for ( k=0; k<-ier; ++k ) {
1637 printf("%" MMG5_PRId " ",MMG3D_indElt(mesh,invalidTets[k]));
1638 }
1639 printf("\n\n");
1640 }
1641 MMG5_SAFE_FREE(invalidTets);
1642 }
1643 disp->npi = disp->np;
1644
1645 if ( (ier > 0) && mesh->info.optim ) {
1646 if ( !MMG3D_doSol(mesh,met) ) {
1647 if ( !MMG5_unscaleMesh(mesh,met,disp) ) _LIBMMG5_RETURN(mesh,met,disp,MMG5_STRONGFAILURE);
1649 }
1650 }
1651
1652 chrono(OFF,&(ctim[3]));
1653 printim(ctim[3].gdif,stim);
1654 if ( mesh->info.imprim > 0 ) {
1655 fprintf(stdout," -- PHASE 2 COMPLETED. %s\n",stim);
1656 }
1657
1658 /* End with a classical remeshing stage, provided mesh->info.lag >= 1 */
1659 if ( (ier > 0) && (mesh->info.lag >= 1) ) {
1660 chrono(ON,&(ctim[4]));
1661 if ( mesh->info.imprim > 0 ) {
1662 fprintf(stdout,"\n -- PHASE 3 : MESH IMPROVEMENT\n");
1663 }
1664
1665 /* renumerotation if available */
1666 if ( !MMG5_scotchCall(mesh,met,disp,NULL) )
1667 {
1668 if ( !MMG5_unscaleMesh(mesh,met,disp) ) _LIBMMG5_RETURN(mesh,met,disp,MMG5_STRONGFAILURE);
1670 }
1671
1672#ifdef MMG_PATTERN
1673 if ( !MMG5_mmg3d1_pattern(mesh,met,NULL) ) {
1674 if ( !(mesh->adja) && !MMG3D_hashTetra(mesh,1) ) {
1675 fprintf(stderr,"\n ## Hashing problem. Invalid mesh.\n");
1677 }
1678 if ( !MMG5_unscaleMesh(mesh,met,disp) ) _LIBMMG5_RETURN(mesh,met,disp,MMG5_STRONGFAILURE);
1680 }
1681#else
1682 if ( !MMG5_mmg3d1_delone(mesh,met,NULL) ) {
1683 if ( !(mesh->adja) && !MMG3D_hashTetra(mesh,1) ) {
1684 fprintf(stderr,"\n ## Hashing problem. Invalid mesh.\n");
1686 }
1687 if ( !MMG5_unscaleMesh(mesh,met,disp) ) _LIBMMG5_RETURN(mesh,met,disp,MMG5_STRONGFAILURE);
1689 }
1690#endif
1691
1692 chrono(OFF,&(ctim[4]));
1693 printim(ctim[4].gdif,stim);
1694 if ( mesh->info.imprim > 0 ) {
1695 fprintf(stdout," -- PHASE 3 COMPLETED. %s\n",stim);
1696 }
1697 }
1698
1699 /* last renum to give back a good numbering to the user */
1700 if ( !MMG5_scotchCall(mesh,met,disp,NULL) )
1701 {
1702 if ( !MMG5_unscaleMesh(mesh,met,disp) ) _LIBMMG5_RETURN(mesh,met,disp,MMG5_STRONGFAILURE);
1704 }
1705
1706 /* save file */
1707 if ( !MMG3D_outqua(mesh,met) ) {
1708 if ( !MMG5_unscaleMesh(mesh,met,disp) ) {
1709 disp->npi = disp->np;
1711 }
1713 }
1714
1715 if ( mesh->info.imprim > 1 && met->m )
1716 MMG3D_prilen(mesh,met,1);
1717
1718 chrono(ON,&(ctim[1]));
1719 if ( mesh->info.imprim > 0 ) fprintf(stdout,"\n -- MESH PACKED UP\n");
1720 if ( !MMG5_unscaleMesh(mesh,met,disp) ) {
1721 disp->npi = disp->np;
1723 }
1724 if ( !MMG3D_packMesh(mesh,met,disp) ) {
1725 disp->npi = disp->np;
1727 }
1728
1729 chrono(OFF,&(ctim[1]));
1730
1731 chrono(OFF,&ctim[0]);
1732 printim(ctim[0].gdif,stim);
1733 if ( mesh->info.imprim >= 0 ) {
1734 fprintf(stdout,"\n MMG3DMOV: ELAPSED TIME %s\n",stim);
1735 fprintf(stdout,"\n %s\n END OF MODULE MMG3D\n %s\n\n",MG_STR,MG_STR);
1736 }
1737 disp->npi = disp->np;
1739}
1740
1746 MMG5_bezierCP = MMG5_mmg3dBezierCP;
1747 MMG5_chkmsh = MMG5_mmg3dChkmsh;
1748 MMG5_indPt = MMG3D_indPt;
1749 MMG5_indElt = MMG3D_indElt;
1750 MMG5_grad2met_ani = MMG5_grad2metSurf;
1751 MMG5_grad2metreq_ani = MMG5_grad2metSurfreq;
1752 MMG5_solTruncature_ani = MMG5_3dSolTruncature_ani;
1753
1754#ifdef USE_SCOTCH
1755 MMG5_renumbering = MMG5_mmg3dRenumbering;
1756#endif
1757}
int ier
MMG5_pMesh MMG5_pSol * sol
if(!ier) exit(EXIT_FAILURE)
MMG5_pMesh * mesh
int MMG3D_analys(MMG5_pMesh mesh)
Definition: analys_3d.c:1324
MMG5_int MMG5_grad2metSurf(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTria pt, MMG5_int np1, MMG5_int np2)
Definition: anisosiz.c:975
int MMG5_grad2metSurfreq(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTria pt, MMG5_int npmaster, MMG5_int npslave)
Definition: anisosiz.c:1951
int MMG5_mmg3dBezierCP(MMG5_pMesh mesh, MMG5_Tria *pt, MMG5_pBezier pb, int8_t ori)
Definition: bezier_3d.c:326
int MMG5_mmg3dChkmsh(MMG5_pMesh mesh, int severe, MMG5_int base)
Definition: chkmsh_3d.c:454
void tminit(mytime *t, int maxtim)
Initialize mytime object.
Definition: chrono.c:120
void printim(double elps, char *stim)
Print real time.
Definition: chrono.c:160
void chrono(int cmode, mytime *ptt)
Function to measure time.
Definition: chrono.c:49
#define TIMEMAX
int MMG3D_hashTetra(MMG5_pMesh mesh, int pack)
Definition: hash_3d.c:122
int MMG5_hNew(MMG5_pMesh mesh, MMG5_HGeom *hash, MMG5_int hsiz, MMG5_int hmax)
Definition: hash_3d.c:994
int MMG3D_hashPrism(MMG5_pMesh mesh)
Definition: hash_3d.c:238
int MMG5_chkBdryTria(MMG5_pMesh mesh)
Definition: hash_3d.c:1439
int MMG5_hEdge(MMG5_pMesh mesh, MMG5_HGeom *hash, MMG5_int a, MMG5_int b, MMG5_int ref, int16_t tag)
Definition: hash_3d.c:951
int MMG3D_pack_sol(MMG5_pMesh mesh, MMG5_pSol sol)
Definition: libmmg3d.c:713
int MMG3D_mark_packedPoints(MMG5_pMesh mesh, MMG5_int *np, MMG5_int *nc)
Definition: libmmg3d.c:541
int MMG3D_packMesh(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_pSol met)
Definition: libmmg3d.c:905
void MMG3D_Free_topoTables(MMG5_pMesh mesh)
Definition: libmmg3d.c:67
void MMG3D_Set_commonFunc(void)
Definition: libmmg3d.c:1745
int MMG3D_pack_tetra(MMG5_pMesh mesh)
Definition: libmmg3d.c:627
int MMG3D_mmg3dls(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_pSol umet)
Definition: libmmg3d.c:1192
int MMG3D_mmg3dlib(MMG5_pMesh mesh, MMG5_pSol met)
Definition: libmmg3d.c:975
int MMG3D_update_eltsVertices(MMG5_pMesh mesh)
Definition: libmmg3d.c:807
MMG5_int MMG3D_bdryBuild(MMG5_pMesh mesh)
Definition: libmmg3d.c:99
MMG5_int MMG3D_pack_points(MMG5_pMesh mesh)
Definition: libmmg3d.c:853
int MMG3D_pack_tetraAndAdja(MMG5_pMesh mesh)
Definition: libmmg3d.c:572
#define MMG5_RETURN_AND_PACK(mesh, met, sol, val)
Definition: libmmg3d.c:51
int MMG3D_mmg3dmov(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol disp)
Definition: libmmg3d.c:1475
int MMG3D_pack_prismsAndQuads(MMG5_pMesh mesh)
Definition: libmmg3d.c:667
int MMG3D_pack_pointArray(MMG5_pMesh mesh)
Definition: libmmg3d.c:749
void MMG3D_unset_reqBoundaries(MMG5_pMesh mesh)
Definition: libmmg3d.c:875
API headers for the mmg3d library.
LIBMMG3D_EXPORT int MMG3D_Clean_isoSurf(MMG5_pMesh mesh)
LIBMMG3D_EXPORT int MMG3D_Set_constantSize(MMG5_pMesh mesh, MMG5_pSol met)
LIBMMG3D_EXPORT int(* MMG3D_doSol)(MMG5_pMesh mesh, MMG5_pSol met)
Definition: mmg3dexterns.c:11
LIBMMG3D_EXPORT void MMG3D_setfunc(MMG5_pMesh mesh, MMG5_pSol met)
int MMG5_mmg3dRenumbering(int, MMG5_pMesh, MMG5_pSol, MMG5_pSol, MMG5_int *)
Definition: librnbg_3d.c:114
int MMG5_mmg3d3(MMG5_pMesh, MMG5_pSol, MMG5_pSol, MMG5_int **)
Definition: mmg3d3.c:583
void MMG3D_keep_only1Subdomain(MMG5_pMesh mesh, int nsd)
Definition: tools_3d.c:1391
void MMG5_freeXPrisms(MMG5_pMesh mesh)
Definition: zaldy_3d.c:378
void MMG3D_delPt(MMG5_pMesh mesh, MMG5_int ip)
Definition: zaldy_3d.c:80
void MMG5_freeXTets(MMG5_pMesh mesh)
Definition: zaldy_3d.c:359
MMG5_int MMG3D_indPt(MMG5_pMesh mesh, MMG5_int kp)
Definition: tools_3d.c:916
int MMG3D_mmg3d2(MMG5_pMesh, MMG5_pSol, MMG5_pSol)
Definition: mmg3d2.c:2208
int MMG3D_delElt(MMG5_pMesh mesh, MMG5_int iel)
Definition: zaldy_3d.c:122
int MMG3D_tetraQual(MMG5_pMesh mesh, MMG5_pSol met, int8_t metRidTyp)
Definition: quality_3d.c:50
static const uint8_t MMG5_iare[6][2]
vertices of extremities of the edges of the tetra
int MMG5_mmg3d1_pattern(MMG5_pMesh, MMG5_pSol, MMG5_int *)
int MMG3D_outqua(MMG5_pMesh mesh, MMG5_pSol met)
Definition: quality_3d.c:764
int MMG5_mmg3d1_delone(MMG5_pMesh, MMG5_pSol, MMG5_int *)
int MMG3D_inqua(MMG5_pMesh mesh, MMG5_pSol met)
Definition: quality_3d.c:631
int MMG3D_prilen(MMG5_pMesh mesh, MMG5_pSol met, int8_t)
Definition: quality_3d.c:341
static void MMG5_warnOrientation(MMG5_pMesh mesh)
MMG5_int MMG3D_indElt(MMG5_pMesh mesh, MMG5_int kel)
Definition: tools_3d.c:860
LIBMMG_CORE_EXPORT int MMG5_unscaleMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol ls)
Definition: scalem.c:689
LIBMMG_CORE_EXPORT int MMG5_scaleMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol ls)
Definition: scalem.c:649
#define MMG5_STRONGFAILURE
Definition: libmmgtypes.h:59
#define MMG5_LOWFAILURE
Definition: libmmgtypes.h:51
#define MMG5_SUCCESS
Definition: libmmgtypes.h:43
int MMG5_scotchCall(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol fields, MMG5_int *permNodGlob)
Definition: librnbg.c:231
void MMG5_version(MMG5_pMesh mesh, char *dim)
Functions needed by libraries API.
Definition: libtools.c:43
#define MG_REQ
#define MG_GEO
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
int MMG5_3dSolTruncature_ani(MMG5_pMesh mesh, MMG5_pSol met)
Definition: scalem.c:453
#define _LIBMMG5_RETURN(mesh, sol, met, val)
#define MG_EOK(pt)
#define MG_NUL
static void MMG5_excfun(int sigid)
#define MMG5_ADD_MEM(mesh, size, message, law)
#define MG_EDG(tag)
#define MG_VOK(ppt)
static void MMG5_warnScotch(MMG5_pMesh mesh)
#define MG_CRN
#define MG_BDY
#define MG_NOSURF
#define MG_NOM
#define MG_STR
#define MG_REF
#define MMG5_SAFE_FREE(ptr)
#define MMG5_DEL_MEM(mesh, ptr)
Structure to store edges of a MMG mesh.
Definition: libmmgtypes.h:305
MMG5_int b
Definition: libmmgtypes.h:306
MMG5_int ref
Definition: libmmgtypes.h:307
int16_t tag
Definition: libmmgtypes.h:310
MMG5_int a
Definition: libmmgtypes.h:306
MMG5_int max
Definition: libmmgtypes.h:576
MMG5_hgeom * geom
Definition: libmmgtypes.h:575
int8_t iso
Definition: libmmgtypes.h:534
int8_t ddebug
Definition: libmmgtypes.h:532
double hsiz
Definition: libmmgtypes.h:518
int8_t isosurf
Definition: libmmgtypes.h:535
int8_t setfem
Definition: libmmgtypes.h:536
uint8_t optimLES
Definition: libmmgtypes.h:546
int8_t lag
Definition: libmmgtypes.h:540
MMG5_int nsd
Definition: libmmgtypes.h:522
uint8_t optim
Definition: libmmgtypes.h:546
int8_t fem
Definition: libmmgtypes.h:539
MMG mesh structure.
Definition: libmmgtypes.h:605
size_t memCur
Definition: libmmgtypes.h:607
MMG5_pQuad quadra
Definition: libmmgtypes.h:648
MMG5_int nc1
Definition: libmmgtypes.h:615
MMG5_Info info
Definition: libmmgtypes.h:651
MMG5_int * adjapr
Definition: libmmgtypes.h:632
MMG5_int ne
Definition: libmmgtypes.h:612
MMG5_pPoint point
Definition: libmmgtypes.h:641
MMG5_int * adja
Definition: libmmgtypes.h:624
MMG5_pPrism prism
Definition: libmmgtypes.h:645
MMG5_int nquad
Definition: libmmgtypes.h:613
MMG5_int nemax
Definition: libmmgtypes.h:612
MMG5_int npmax
Definition: libmmgtypes.h:612
MMG5_pxPoint xpoint
Definition: libmmgtypes.h:642
MMG5_int nenil
Definition: libmmgtypes.h:622
MMG5_HGeom htab
Definition: libmmgtypes.h:650
MMG5_int xp
Definition: libmmgtypes.h:620
MMG5_int nei
Definition: libmmgtypes.h:612
MMG5_pTetra tetra
Definition: libmmgtypes.h:643
MMG5_int nt
Definition: libmmgtypes.h:612
MMG5_pTria tria
Definition: libmmgtypes.h:647
MMG5_int np
Definition: libmmgtypes.h:612
MMG5_pEdge edge
Definition: libmmgtypes.h:649
MMG5_int nti
Definition: libmmgtypes.h:612
MMG5_int npi
Definition: libmmgtypes.h:612
MMG5_pxTetra xtetra
Definition: libmmgtypes.h:644
MMG5_int nai
Definition: libmmgtypes.h:612
MMG5_int nprism
Definition: libmmgtypes.h:613
MMG5_int na
Definition: libmmgtypes.h:612
MMG5_int npnil
Definition: libmmgtypes.h:621
Structure to store points of a MMG mesh.
Definition: libmmgtypes.h:270
double n[3]
Definition: libmmgtypes.h:272
int16_t tag
Definition: libmmgtypes.h:284
MMG5_int tmp
Definition: libmmgtypes.h:280
MMG5_int xp
Definition: libmmgtypes.h:279
MMG5_int ref
Definition: libmmgtypes.h:278
MMG5_int v[6]
Definition: libmmgtypes.h:463
MMG5_int v[4]
Definition: libmmgtypes.h:367
MMG5_int npi
Definition: libmmgtypes.h:667
double * m
Definition: libmmgtypes.h:671
MMG5_int np
Definition: libmmgtypes.h:665
MMG5_int v[4]
Definition: libmmgtypes.h:403
MMG5_int xt
Definition: libmmgtypes.h:407
Cell of the hash table of geom edges.
Definition: libmmgtypes.h:562
MMG5_int ref
Definition: libmmgtypes.h:565
int16_t tag
Definition: libmmgtypes.h:567
MMG5_int b
Definition: libmmgtypes.h:564
MMG5_int a
Definition: libmmgtypes.h:563
double n1[3]
Definition: libmmgtypes.h:295
int16_t tag[6]
Definition: libmmgtypes.h:425
MMG5_int edg[6]
Definition: libmmgtypes.h:421
Chrono object.