Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
colver_3d.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
38
39extern int8_t ddb;
40
43int MMG5_chkcol_int(MMG5_pMesh mesh,MMG5_pSol met,MMG5_int k,int8_t iface,
44 int8_t iedg,int64_t *list,int ilist,int8_t typchk) {
45 MMG5_pTetra pt,pt0;
46 MMG5_pPoint p0;
47 double calold,calnew,caltmp,ll,lon;
48 MMG5_int j,iel,nq,nr;
49 int8_t i,jj,ip,iq;
50
51 iq = MMG5_idir[iface][MMG5_iprv2[iedg]];
52 pt = &mesh->tetra[k];
53 pt0 = &mesh->tetra[0];
54 nq = pt->v[iq];
55
56 lon = 1.6;
57 if ( typchk == 2 && met->m ) {
58 lon = MMG5_lenedg(mesh,met,MMG5_iarf[iface][iedg],pt);
59
60 if ( !lon ) return 0;
61 /*on cherche a se rapprocher de 1*/
62 //lon = MG_MAX(0.8/lon,1.6);// test ok but less good than the next one
63 lon = MG_MAX(2.-lon,1.6);
64 }
65
66 calold = calnew = DBL_MAX;
67 for (j=0; j<ilist; j++) {
68 iel = list[j] / 4;
69 ip = list[j] % 4;
70 pt = &mesh->tetra[iel];
71
72 /* exclude elements from shell */
73 for (jj=0; jj<4; jj++) if ( pt->v[jj] == nq ) break;
74 if ( jj < 4 ) continue;
75 memcpy(pt0,pt,sizeof(MMG5_Tetra));
76 /* Update edges tag for pt0 (needed by lenedg). */
77
78 /* prevent from recreating internal edge between boundaries */
79 if ( mesh->info.fem==typchk ) {
80 p0 = &mesh->point[nq];
81 if ( (p0->tag & MG_BDY) && !(p0->tag & MG_PARBDY) ) {
82 i = ip;
83 for (jj=0; jj<3; jj++) {
84 i = MMG5_inxt3[i];
85 p0 = &mesh->point[pt->v[i]];
86 if ( (p0->tag & MG_BDY) && !(p0->tag & MG_PARBDY) ){
87 return 0;
88 }
89 }
90
91 /* Prevent from creating a tetra with 4 bdy vertices */
92 // Algiane (2022) this test is useless I think (because we forbid the
93 // creation of internal edges between boundary points)
94#ifndef NDEBUG
95 i = ip;
96 nr = 0;
97 for (jj=0; jj<3; jj++) {
98 i = MMG5_inxt3[i];
99 p0 = &mesh->point[pt->v[i]];
100 if ( (p0->tag & MG_BDY) && !(p0->tag & MG_PARBDY) ) ++nr;
101 }
102 if ( nr==3 ) {
103 assert ( 0 && "Uncomment this test, it is not useless");
104 return 0;
105 }
106#endif
107 }
108 }
109 else {
110 /* In aniso : prevent from creating a tetra with 4 ridges vertices or
111 * internal edges between two ridges (unable to split it after because of
112 * the ridge metric) */
113 if ( met->size==6 ) {
114 p0 = &mesh->point[nq];
115
116 if ( (p0->tag & MG_GEO) && !MG_SIN_OR_NOM(p0->tag) ) {
117 i = ip;
118 for (jj=0; jj<3; jj++) {
119 i = MMG5_inxt3[i];
120 p0 = &mesh->point[pt->v[i]];
121 if ( p0->tag & MG_GEO && !MG_SIN_OR_NOM(p0->tag) ) {
122 return 0;
123 }
124 }
125
126 // Algiane (2022) this test is useless because we forbid the creation of
127 // internal edges between boundary points but we can keep it in case we
128 // comment the previous test
129#ifndef NDEBUG
130 i = ip;
131 nr = 0;
132 for (jj=0; jj<3; jj++) {
133 i = MMG5_inxt3[i];
134 p0 = &mesh->point[pt->v[i]];
135 if ( p0->tag & MG_GEO && !MG_SIN_OR_NOM(p0->tag) ) {
136 ++nr;
137 }
138 }
139 if ( nr==3 ) {
140 assert ( 0 && "Uncomment this test, it is not useless");
141 return 0;
142 }
143#endif
144 }
145 }
146 }
147
148 pt0->v[ip] = nq;
149
150 calold = MG_MIN(calold,pt->qual);
151 if ( typchk==1 && met->m && met->size > 1 )
152 caltmp = MMG5_caltet33_ani(mesh,met,pt0);
153 else
154 caltmp = MMG5_orcal(mesh,met,0);
155
156 if ( caltmp < MMG5_NULKAL ) return 0;
157 calnew = MG_MIN(calnew,caltmp);
158 /* check length */
159 if ( typchk == 2 && met->m ) {
160 for (jj=0; jj<6; jj++) {
161 /* Rough evaluation of edge length (doesn't take into account if some of
162 * the modified edges of pt0 are boundaries): for a more precise
163 * computation, we need to update the edge tags of pt0. */
164 ll = MMG5_lenedgspl(mesh,met,jj,pt0);
165 if ( (!ll) || (ll > lon) )//LOPTL too small, we need to put greater than 1.41
166 return 0;
167 }
168 }
169 }
170 if ( calold < MMG5_EPSOK && calnew <= calold ) return 0;
171 else if ( calnew < MMG5_EPSOK || calnew < 0.3*calold ) return 0;
172
173 return ilist;
174}
175
192static inline
193MMG5_int MMG3D_unfold_shell(MMG5_pMesh mesh,MMG5_int start,MMG5_int end, MMG5_int na, MMG5_int nb,MMG5_int piv,
194 MMG5_int *iel,int8_t *iopp) {
195 MMG5_pTetra pt;
196 MMG5_int adj,*adja;
197 int8_t i,ipiv,isface;
198
199 adj = start;
200 do {
201 *iel = adj;
202 pt = &mesh->tetra[*iel];
203 adja = &mesh->adja[4*(*iel-1)+1];
204
205 /* Identification of edge number in tetra iel */
206 if ( !MMG3D_findEdge(mesh,pt,*iel,na,nb,1,NULL,&i) ) return -1;
207
208 /* set sense of travel */
209 if ( pt->v[ MMG5_ifar[i][0] ] == piv ) {
210 adj = adja[ MMG5_ifar[i][0] ] / 4;
211 ipiv = MMG5_ifar[i][1];
212 *iopp = MMG5_ifar[i][0];
213 piv = pt->v[ipiv];
214 }
215 else {
216 adj = adja[ MMG5_ifar[i][1] ] / 4;
217 ipiv = MMG5_ifar[i][0];
218 *iopp = MMG5_ifar[i][1];
219 piv = pt->v[ipiv];
220 }
221
222 isface = 0;
223 if ( pt->xt ) {
224 isface = (MG_BDY & mesh->xtetra[pt->xt].ftag[*iopp]);
225 }
226 }
227 while ( adj && ( adj != end ) && !isface );
228
229 return piv;
230}
231
250static int
251MMG5_topchkcol_bdy(MMG5_pMesh mesh,MMG5_int k,int iface,int8_t iedg,MMG5_int *lists,
252 int ilists) {
253 MMG5_pTetra pt;
254 MMG5_pxTetra pxt;
255 double n0[3],n1[3],devnew;
256 MMG5_int nump,numq,iel,jel,jel1,nap,nbp,naq,nbq,nro;
257 int8_t ip,iq,iopp,i,j,j1,jface,jface1;
258
259 pt = &mesh->tetra[k];
260 ip = MMG5_idir[iface][MMG5_inxt2[iedg]];
261 iq = MMG5_idir[iface][MMG5_iprv2[iedg]];
262 nump = pt->v[ip];
263 numq = pt->v[iq];
264
265 /* Surface ball has been enumerated as f1,...,f2 - f1,f2 = both triangles of
266 * surface shell, point nap, facing the first vanishing face in surface ball
267 * of p */
268 nro = pt->v[MMG5_idir[iface][iedg]];
269
270 jel = lists[1] / 4;
271 jface = lists[1] % 4;
272
273 pt = &mesh->tetra[jel];
274 for (j=0; j<3; j++) {
275 i = MMG5_idir[jface][j];
276 if ( pt->v[i] != nump && pt->v[i] != nro ) break;
277 }
278 assert(j<3);
279
280 nap = pt->v[i];
281
282 /* Unfold shell of (numq,nro), starting from (jel,jface), with pivot nump */
283 naq = MMG3D_unfold_shell(mesh,k,k,numq,nro,nump,&iel,&iopp);
284
285 if ( naq < 0 ) {
286 return -1;
287 }
288 else if ( nap == naq ) {
289 return 0;
290 }
291
292 assert ( mesh->tetra[k].xt && "initial tetra is not boundary");
293 pxt = &mesh->xtetra[mesh->tetra[k].xt];
294
295 if ( !( MG_GEO_OR_NOM(pxt->tag[MMG5_iarf[iface][MMG5_inxt2[iedg]]]) ||
296 MG_GEO_OR_NOM(pxt->tag[MMG5_iarf[iface][MMG5_iprv2[iedg]]]) ) ) {
297
298 /* Check the normal deviation between the boundary faces sharing the edge
299 * numq (or nump)-nro */
300 if ( !MMG5_norpts (mesh,numq,nro,nap,n0) ) return 0;
301 if ( !MMG5_norface(mesh,iel,iopp ,n1) ) return 0;
302
303 devnew = n0[0]*n1[0] + n0[1]*n1[1] + n0[2]*n1[2];
304 if ( devnew < mesh->info.dhd ) return 0;
305 }
306
307 /* Point nbp, facing the second vanishing face in surface ball of p */
308 jel1 = lists[ilists-1] / 4;
309 jface1 = lists[ilists-1] % 4;
310 pt = &mesh->tetra[jel1];
311 for (j1=0; j1<3; j1++) {
312 i = MMG5_idir[jface1][j1];
313 if ( pt->v[i] != nump && pt->v[i] != numq ) break;
314 }
315 assert(j1<3);
316
317 nro = pt->v[i];
318 jel = lists[ilists-2] / 4;
319 jface = lists[ilists-2] % 4;
320 pt = &mesh->tetra[jel];
321 for (j=0; j<3; j++) {
322 i = MMG5_idir[jface][j];
323 if ( pt->v[i] != nump && pt->v[i] != nro ) break;
324 }
325 assert(j<3);
326
327 nbp = pt->v[i];
328
329 /* Unfold shell of (numq,nro), starting from (jel,jface), with pivot nump */
330 nbq = MMG3D_unfold_shell(mesh,lists[ilists-1]/4,k,numq,nro,nump,&iel,&iopp);
331
332 if ( nbq < 0 ) {
333 return -1;
334 }
335 else if ( nbp == nbq ) {
336 return 0;
337 }
338
339 pxt = &mesh->xtetra[mesh->tetra[jel1].xt];
340
341 if ( !( MG_GEO_OR_NOM(pxt->tag[MMG5_iarf[jface1][MMG5_iprv2[j1]]]) ||
342 MG_GEO_OR_NOM(pxt->tag[MMG5_iarf[jface1][MMG5_inxt2[j1]]]) ) ) {
343
344 /* Check the normal deviation between the boundary faces sharing the edge
345 * numq (or nump)-nro */
346 if ( !MMG5_norpts (mesh,nro,numq,nbp,n0) ) return 0;
347 if ( !MMG5_norface(mesh,iel,iopp ,n1) ) return 0;
348
349 devnew = n0[0]*n1[0] + n0[1]*n1[1] + n0[2]*n1[2];
350 if ( devnew < mesh->info.dhd ) return 0;
351 }
352
353 return 1;
354}
355
376static inline
377int MMG3D_get_shellEdgeTag_oneDir(MMG5_pMesh mesh,MMG5_int start, MMG5_int na, MMG5_int nb,
378 int16_t *tag,MMG5_int *ref, MMG5_int piv,MMG5_int adj,
379 int8_t *filled) {
380 MMG5_pTetra pt;
381 MMG5_pxTetra pxt;
382 MMG5_int *adja;
383 int8_t i;
384
385 *filled = 0;
386 while ( adj && (adj != start) ) {
387 pt = &mesh->tetra[adj];
388
389 /* identification of edge number in tetra adj */
390 if ( !MMG3D_findEdge(mesh,pt,adj,na,nb,1,NULL,&i) ) return -1;
391
392 /* update edge ref and tag */
393 if ( pt->xt ) {
394 pxt = &mesh->xtetra[pt->xt];
395 *ref = pxt->edg[i];
396 if ( pxt->tag[i] & MG_BDY ) {
397 *tag |= pxt->tag[i];
398 *filled = 1;
399 return adj;
400 }
401 }
402
403 /* set new triangle for travel */
404 adja = &mesh->adja[4*(adj-1)+1];
405 if ( pt->v[ MMG5_ifar[i][0] ] == piv ) {
406 adj = adja[ MMG5_ifar[i][0] ] / 4;
407 piv = pt->v[ MMG5_ifar[i][1] ];
408 }
409 else {
410 adj = adja[ MMG5_ifar[i][1] ] /4;
411 piv = pt->v[ MMG5_ifar[i][0] ];
412 }
413 }
414
415 return adj;
416}
417
431static inline
432int MMG3D_get_shellEdgeTag(MMG5_pMesh mesh,MMG5_int start, int8_t ia,int16_t *tag,MMG5_int *ref) {
433 MMG5_pTetra pt;
434 MMG5_pxTetra pxt;
435 MMG5_int piv,na,nb,adj,*adja;
436 int8_t filled;
437
438 pt = &mesh->tetra[start];
439
440 assert( start >= 1 && MG_EOK(pt) );
441
442 pxt = NULL;
443 na = pt->v[MMG5_iare[ia][0]];
444 nb = pt->v[MMG5_iare[ia][1]];
445
446 if ( pt->xt ) {
447 pxt = &mesh->xtetra[pt->xt];
448 if ( pxt->tag[ia] & MG_BDY ) {
449 *tag |= pxt->tag[ia];
450 *ref = pxt->edg[ia];
451 return 1;
452 }
453 }
454
455 /* Travel in one direction */
456 adja = &mesh->adja[4*(start-1)+1];
457 adj = adja[MMG5_ifar[ia][0]] / 4;
458 piv = pt->v[MMG5_ifar[ia][1]];
459
460 adj = MMG3D_get_shellEdgeTag_oneDir(mesh,start,na,nb,tag,ref,piv,adj,&filled);
461
462 /* If an xtetra has been found or if all shell has been travelled, stop, else,
463 * travel it the other sense */
464 if ( adj > 0 ) {
465 assert ( (adj == start) || filled );
466 return 1;
467 }
468 else if ( adj < 0 ) return 0;
469
470 assert(!adj);
471
472 pt = &mesh->tetra[start];
473 adja = &mesh->adja[4*(start-1)+1];
474 adj = adja[MMG5_ifar[ia][1]] / 4;
475 piv = pt->v[MMG5_ifar[ia][0]];
476
477 adj = MMG3D_get_shellEdgeTag_oneDir(mesh,start,na,nb,tag,ref,piv,adj,&filled);
478
479 if ( adj < 0 ) return 0;
480
481 return 1;
482}
483
515int MMG5_chkcol_bdy(MMG5_pMesh mesh,MMG5_pSol met,MMG5_int k,int8_t iface,
516 int8_t iedg,int64_t *listv,int ilistv,MMG5_int *lists,int ilists,
517 MMG5_int refmin,MMG5_int refplus, int8_t typchk,int isnm,int8_t isnmint) {
518 MMG5_pTetra pt,pt0,pt1;
519 MMG5_pxTetra pxt;
520 MMG5_Tria tt;
521 MMG5_pPar par;
522 double calold,calnew,caltmp,nadja[3],nprvold[3],nprvnew[3],ncurold[3],ncurnew[3];
523 double ps,devold,devnew,hmax,hausd;
524 MMG5_int nump,numq,ndepmin,ndepplus,l,kk,iel;
525 int nr,nbbdy,isloc,iedgeOpp,ipp;
526 int16_t tag;
527 int8_t iopp,iopp2,ia,ip,i,iq,i0,i1,ier,isminp,isplp;
528#ifndef NDEBUG
529 MMG5_pPoint p0;
530#endif
531
532 pt = &mesh->tetra[k];
533 pxt = 0;
534 pt0 = &mesh->tetra[0];
535 ip = MMG5_idir[iface][MMG5_inxt2[iedg]];
536 nump = pt->v[ip];
537 numq = pt->v[MMG5_idir[iface][MMG5_iprv2[iedg]]];
538
539#ifndef NDEBUG
540 p0 = &mesh->point[nump];
541 assert(p0->tag & MG_BDY);
542 assert(p0->xp);
543
544 /* Remove maybe-uninitialized value warning */
545 nprvold[0] = nprvold[1] = nprvold[2] = 0.;
546 ncurold[0] = ncurold[1] = ncurold[2] = 0.;
547 nprvnew[0] = nprvnew[1] = nprvnew[2] = 0.;
548 ncurnew[0] = ncurnew[1] = ncurnew[2] = 0.;
549#endif
550
551 ndepmin = ndepplus = 0;
552 isminp = isplp = 0;
553
554 if ( !isnm ) {
555 refmin = MG_MINUS;
556 refplus = MG_PLUS;
557 }
558
559 if ( !isnmint ) {
560 /* prevent collapse in case surface ball has 3 triangles */
561 if ( ilists <= 2 ) return 0; // ATTENTION, Normalement, avec 2 c est bon !
562
563 /* Surfacic ball is enumerated with first tet having (pq) as edge n° MMG5_iprv2[ip] on face iopp */
564 MMG5_startedgsurfball(mesh,nump,numq,lists,ilists);
565 }
566
567 /* check tetra quality */
568 calold = calnew = DBL_MAX;
569 for (l=0; l<ilistv; l++) {
570 iel = listv[l] / 4;
571 ipp = listv[l] % 4;
572 assert(iel && ipp >=0 && ipp < 4 && "unexpected tetra or local vertex idx");
573
574 pt = &mesh->tetra[iel];
575
576 if ( !isnmint ) {
577 if ( pt->ref == refmin ) isminp = 1;
578 else if ( pt->ref == refplus ) isplp = 1;
579 }
580
581 /* Topological test for tetras of the shell */
582 for (iq=0; iq<4; iq++)
583 if ( pt->v[iq] == numq ) break;
584
585 if ( iq < 4 ) {
586 nbbdy = 0;
587 if ( pt->xt ) {
588 pxt = &mesh->xtetra[pt->xt];
589 }
590
591 for (i=0; i<4; i++) {
592 if ( pt->xt && (pxt->ftag[i] & MG_BDY) ) nbbdy++;
593 }
594
595 if ( nbbdy == 4 ) {
596 /* Only one element in the domain: we don't want to delete it */
597 return 0;
598 }
599 else if ( nbbdy == 3 ) {
600 /* Topological problem triggered when one of the two faces of collapsed
601 edge is the only internal one : closing a part of the domain */
602
603 /* Identification of edge number in tetra iel */
604 if ( !MMG3D_findEdge(mesh,pt,iel,numq,nump,1,NULL,&ia) ) return -1;
605
606 i0 = MMG5_ifar[ia][0];
607 i1 = MMG5_ifar[ia][1];
608 if ( pt->xt && (!(pxt->ftag[i0] & MG_BDY) || !(pxt->ftag[i1] & MG_BDY)) )
609 return 0;
610 }
611
612 /* Now check that the 2 faces identified by collapse are not boundary */
613 if ( pt->xt && (pxt->ftag[ipp] & MG_BDY) && (pxt->ftag[iq] & MG_BDY) )
614 return 0;
615
616 if ( pt->xt ) {
617 for (i=0; i<4; i++) {
618 if ( i==ipp || i==iq ) {
619 continue;
620 }
621
622 /* Avoid surface crimping: check that the collapse doesn't merge 3
623 * bdy edge along a non bdy face: we have to check the edge of each
624 * shell because some MG_BDY tags may be missings due to the creation
625 * of an xtetra during a previous collapse */
626 if ( !(pxt->ftag[i] & MG_BDY) ) {
627 int16_t tag0,tag1,tag2;
628 int ref0,ref1,ref2;
629
630 tag0 = tag1 = tag2 = 0;
631 ref0 = ref1 = ref2 = 0;
632
633 if ( !MMG3D_get_shellEdgeTag(mesh,iel,MMG5_iarf[i][0],&tag0,&ref0) ) {
634 fprintf(stderr,"\n ## Error: %s: 0. unable to get edge info.\n",__func__);
635 return 0;
636 }
637 if ( !MMG3D_get_shellEdgeTag(mesh,iel,MMG5_iarf[i][1],&tag1,&ref1) ) {
638 fprintf(stderr,"\n ## Error: %s: 1. unable to get edge info.\n",__func__);
639 return 0;
640 }
641 if ( !MMG3D_get_shellEdgeTag(mesh,iel,MMG5_iarf[i][2],&tag2,&ref2) ) {
642 fprintf(stderr,"\n ## Error: %s: 2. unable to get edge info.\n",__func__);
643 return 0;
644 }
645 if ( (tag0 & MG_BDY) && (tag1 & MG_BDY) && (tag2 & MG_BDY ) ) {
646 return 0;
647 }
648 }
649 }
650 }
651
652 continue;
653 }
654
655 if ( !isnmint ) {
656 if ( isnm || mesh->info.iso ) {
657 /* Volume test for tetras outside the shell */
658 if ( (!ndepmin) && (pt->ref == refmin) ) {
659 ndepmin = iel;
660 }
661 else if ( (!ndepplus) && (pt->ref == refplus) ) {
662 ndepplus = iel;
663 }
664 }
665 }
666
667 /* Prevent from creating a tetra with 4 ridges metrics in aniso mode */
668 if ( met && met->m && met->size == 6 ) {
669 if ( (mesh->point[numq].tag & MG_GEO) && !MG_SIN_OR_NOM(mesh->point[numq].tag) ) {
670 i = ipp;
671 nr = 0;
672 for (iq=0; iq<3; iq++) {
673 i = MMG5_inxt3[i];
674 if ( (mesh->point[pt->v[i]].tag & MG_GEO) &&
675 !MG_SIN_OR_NOM(mesh->point[pt->v[i]].tag) ) {
676 ++nr;
677 }
678 }
679 if ( nr==3 ) return 0;
680 }
681 }
682
683 memcpy(pt0,pt,sizeof(MMG5_Tetra));
684 pt0->v[ipp] = numq;
685
686 calold = MG_MIN(calold, pt->qual);
687 if ( typchk==1 && met->m && met->size > 1 )
688 caltmp = MMG5_caltet33_ani(mesh,met,pt0);
689 else
690 caltmp = MMG5_orcal(mesh,met,0);
691
692 if ( caltmp < MMG5_NULKAL ) return 0;
693 calnew = MG_MIN(calnew,caltmp);
694 }
695 if ( calold < MMG5_EPSOK && calnew <= calold ) return 0;
696 else if ( calnew < MMG5_EPSOK || calnew < 0.3*calold ) return 0;
697
698 if ( isnmint ) {
699 return ilistv;
700 }
701
702 /* analyze surfacic ball of p */
703 for (l=1; l<ilists-1; l++) {
704 iel = lists[l] / 4;
705 iopp = lists[l] % 4;
706 assert(iel && iopp >=0 && iopp < 4 && "unexpected tetra or local vertex idx");
707
708 pt = &mesh->tetra[iel];
709 pxt = &mesh->xtetra[pt->xt];
710 assert(pt->xt);
711
712 /* retrieve vertex in tetra */
713 for (ip=0; ip<4; ip++)
714 if ( pt->v[ip] == nump ) break;
715
716 assert(ip<4);
717 if ( ip==4 ) return 0;
718
719 memcpy(pt0,pt,sizeof(MMG5_Tetra));
720 pt0->v[ip] = numq;
721
722 if ( !MMG5_norface(mesh,iel,iopp,ncurold) ) return 0;
723 if ( !MMG5_norface(mesh,0,iopp,ncurnew) ) return 0;
724
725 /* check normal flipping */
726 ps = ncurold[0]*ncurnew[0] + ncurold[1]*ncurnew[1] + ncurold[2]*ncurnew[2];
727 if ( ps < 0.0 ) return 0;
728
729 /* check normal deviation between l and its neighbour through the edge ia */
730 ia = MMG5_idirinv[iopp][ip]; /* index of p in tria iopp */
731
732 iedgeOpp = MMG5_iarf[iopp][ia];
733
734 if ( ! ( MG_GEO_OR_NOM(mesh->xtetra[pt->xt].tag[iedgeOpp])) ) {
735
736 ier = MMG3D_normalAdjaTri(mesh,iel,iopp,ia,nadja);
737 if ( ier < 0 ) return -1;
738 else if (!ier ) return 0;
739
740 devnew = nadja[0]*ncurnew[0] + nadja[1]*ncurnew[1] + nadja[2]*ncurnew[2];
741
742 if ( devnew < mesh->info.dhd ) return 0;
743 }
744
745 if ( l == 1 ) {
746 /* check normal deviation between the new first tri of the surfacic ball
747 * and the surfacic triangle facing ip in the old first tri of the
748 * surfacic ball (that vanishes due to the collapse)
749 */
750 kk = lists[0] / 4;
751 iopp2 = lists[0] % 4;
752 pt1 = &mesh->tetra[kk];
753 assert(pt1->xt);
754
755 /* retrieve vertex ip in tria iopp2 */
756 for (ipp=0; ipp<3; ipp++)
757 if ( pt1->v[MMG5_idir[iopp2][ipp]] == nump ) break;
758 assert(ipp<3);
759
760 iedgeOpp = MMG5_iarf[iopp2][ipp];
761
762 if ( ! ( MG_GEO_OR_NOM(mesh->xtetra[pt1->xt].tag[iedgeOpp] ) ) ) {
763 ier = MMG3D_normalAdjaTri(mesh,kk,iopp2,ipp,nadja);
764
765 if ( ier < 0 ) return -1;
766 else if ( !ier ) return 0;
767
768 devnew = nadja[0]*ncurnew[0] + nadja[1]*ncurnew[1] + nadja[2]*ncurnew[2];
769 if ( devnew < mesh->info.dhd ) {
770 return 0;
771 }
772 }
773 }
774 else {
775 /* check normal deviation between l and l-1 */
776 ia = MMG5_iprv2[ia]; /* edge between l-1 and l, in local num of tria */
777 ia = MMG5_iarf[iopp][ia]; /* edge between l-1 and l in local num of tetra */
778
779 if ( !MG_GEO_OR_NOM(mesh->xtetra[pt->xt].tag[ia]) ) {
780 devold = nprvold[0]*ncurold[0] + nprvold[1]*ncurold[1] + nprvold[2]*ncurold[2];
781 devnew = nprvnew[0]*ncurnew[0] + nprvnew[1]*ncurnew[1] + nprvnew[2]*ncurnew[2];
782
783 if ( devold < MMG5_ANGEDG ) {
784 if ( devnew < devold ) {
785 return 0;
786 }
787 }
788 else if ( devnew < MMG5_ANGEDG ) {
789 return 0;
790 }
791 }
792 }
793
794 /* check Hausdorff distance to geometric support */
795 MMG5_tet2tri(mesh,iel,iopp,&tt);
796 if ( l == 1 ) {
797 for (i=0; i<3; i++) {
798 if ( tt.v[i] == nump ) break;
799 }
800
801 assert(i<3);
802 if ( i==3 ) return 0;
803
804 /* Index of the third point of the first collapsed triangle */
805 i = MMG5_inxt2[i];
806 ia = MMG5_inxt2[i];
807 tag = pxt->tag[MMG5_iarf[iopp][ia]];
808 tt.tag[ia] = MG_MAX(tt.tag[ia],tag);
809 }
810 else if ( l == ilists-2 ) {
811 for (i=0; i<3; i++) {
812 if ( tt.v[i] == nump ) break;
813 }
814 assert(i<3);
815 /* Index of the third point of the second collapsed triangle */
816 i = MMG5_iprv2[i];
817 ia = MMG5_iprv2[i];
818 tag = pxt->tag[MMG5_iarf[iopp][ia]];
819 tt.tag[ia] = MG_MAX(tt.tag[ia],tag);
820 }
821
822 for (i=0; i<3; i++) {
823 if ( tt.v[i] == nump ) break;
824 }
825 assert(i<3);
826 if ( i==3 ) return 0;
827
828 tt.v[i] = numq;
829
830 /* Local parameters for tt and iel */
831 hmax = mesh->info.hmax;
832 hausd = mesh->info.hausd;
833 isloc = 0;
834 if ( mesh->info.parTyp & MG_Tetra ) {
835 for ( kk=0; kk<mesh->info.npar; ++kk ) {
836 par = &mesh->info.par[kk];
837
838 if ( par->elt != MMG5_Tetrahedron ) continue;
839 if ( par->ref != pt->ref ) continue;
840
841 hmax = par->hmax;
842 hausd = par->hausd;
843 isloc = 1;
844 break;
845 }
846 }
847 if ( mesh->info.parTyp & MG_Tria ) {
848 if ( isloc ) {
849 for ( kk=0; kk<mesh->info.npar; ++kk ) {
850 par = &mesh->info.par[kk];
851
852 if ( par->elt != MMG5_Triangle ) continue;
853 if ( par->ref != tt.ref ) continue;
854
855 hmax = MG_MIN(hmax,par->hmax);
856 hausd = MG_MIN(hausd,par->hausd);
857 break;
858 }
859 }
860 else {
861 for ( kk=0; kk<mesh->info.npar; ++kk ) {
862 par = &mesh->info.par[kk];
863
864 if ( par->elt != MMG5_Triangle ) continue;
865 if ( par->ref != tt.ref ) continue;
866
867 hmax = par->hmax;
868 hausd = par->hausd;
869 isloc = 1;
870 break;
871 }
872 }
873 }
874
875 /* If some edges must be splitted, refuse the collapse */
876 if ( MMG5_chkedg(mesh,&tt,MG_GET(pxt->ori,iopp),hmax,hausd,isloc) > 0 ) return 0;
877
878 memcpy(nprvold,ncurold,3*sizeof(double));
879 memcpy(nprvnew,ncurnew,3*sizeof(double));
880 }
881
882 /* check normal deviation between the new last tri of the surfacic ball
883 * and the surfacic triangle facing ip in the old last tri of the
884 * surfacic ball (that vanishes due to the collapse)
885 */
886 kk = lists[ilists-1] / 4;
887 iopp2 = lists[ilists-1] % 4;
888 pt1 = &mesh->tetra[kk];
889 assert(pt1->xt);
890
891 /* retrieve vertex ip in tria iopp2 */
892 for (ipp=0; ipp<3; ipp++)
893 if ( pt1->v[MMG5_idir[iopp2][ipp]] == nump ) break;
894 assert(ipp<3);
895
896 iedgeOpp = MMG5_iarf[iopp2][ipp];
897
898 if ( ! ( MG_GEO_OR_NOM(mesh->xtetra[pt1->xt].tag[iedgeOpp]) ) ) {
899
900 ier = MMG3D_normalAdjaTri(mesh,kk,iopp2,ipp,nadja);
901 if ( ier < 0 ) return -1;
902 else if ( !ier ) return 0;
903
904 devnew = nadja[0]*ncurnew[0] + nadja[1]*ncurnew[1] + nadja[2]*ncurnew[2];
905
906 if ( devnew < mesh->info.dhd ) {
907 return 0;
908 }
909 }
910
911 if ( isnm || mesh->info.iso ) {
912 ier = MMG3D_chkmanicoll(mesh,k,iface,iedg,ndepmin,ndepplus,refmin,refplus,isminp,isplp);
913 if ( !ier ) return 0;
914 }
915 else {
916 /* Topological check for surface ball */
917 ier = MMG5_topchkcol_bdy(mesh,k,iface,iedg,lists,ilists);
918 if ( ier<0 ) return -1;
919 else if ( !ier ) return 0;
920 }
921
922 return ilistv;
923}
924
939static inline
940void MMG3D_update_edgeTag(MMG5_pTetra pt,MMG5_pxTetra pxt,MMG5_int np, MMG5_int nq,
941 uint8_t ip, MMG5_pTetra pt1,MMG5_pxTetra pxt1,
942 uint8_t voyp) {
943
944 int i,j;
945 MMG5_int p0,p1;
946 uint8_t ia,iav;
947 int16_t tag,tag1;
948
949 /* update tags for edges */
950 for ( j=0; j<3; j++ ) {
951 ia = MMG5_iarf[ip][j];
952 p0 = pt->v[MMG5_iare[ia][0]];
953 p1 = pt->v[MMG5_iare[ia][1]];
954 if ( pxt->tag[ia] ) {
955 for ( i=0; i<3; i++ ) {
956 iav=MMG5_iarf[voyp][i];
957 if ( p0==nq ) {
958 if ( ((pt1->v[MMG5_iare[iav][0]]==np) && (pt1->v[MMG5_iare[iav][1]]==p1)) ||
959 ((pt1->v[MMG5_iare[iav][0]]==p1) && (pt1->v[MMG5_iare[iav][1]]==np)) )
960 break;
961 }
962 else if ( p1==nq ) {
963 if ( ((pt1->v[MMG5_iare[iav][0]]==np ) && (pt1->v[MMG5_iare[iav][1]]==p0)) ||
964 ((pt1->v[MMG5_iare[iav][0]]==p0) && (pt1->v[MMG5_iare[iav][1]]==np )) )
965 break;
966 }
967 else {
968 if ( ((pt1->v[MMG5_iare[iav][0]]==p0) && (pt1->v[MMG5_iare[iav][1]]==p1)) ||
969 ((pt1->v[MMG5_iare[iav][0]]==p1) && (pt1->v[MMG5_iare[iav][1]]==p0)) )
970 break;
971 }
972 }
973 assert(i!=3);
974 tag1 = pxt1->tag[iav];
975 tag = pxt->tag[ia];
976 pxt1->tag[iav] |= pxt->tag[ia];
977 if( ((tag1 & MG_REQ) && !(tag1 & MG_NOSURF)) ||
978 (( tag & MG_REQ) && !( tag & MG_NOSURF)) )
979 pxt1->tag[iav] &= ~MG_NOSURF;
980 pxt1->edg[iav] = MG_MAX ( pxt1->edg[iav],pxt->edg[ia] );
981 }
982 }
983}
984
1001static inline
1002MMG5_int MMG3D_update_shellEdgeTag_oneDir(MMG5_pMesh mesh,MMG5_int start, MMG5_int na, MMG5_int nb,
1003 int16_t tag,MMG5_int ref, MMG5_int piv,MMG5_int adj) {
1004 MMG5_pTetra pt;
1005 MMG5_pxTetra pxt;
1006 MMG5_int *adja;
1007 int16_t xtag;
1008 int8_t i;
1009
1010 assert ( tag & MG_BDY && "Unexpected non boundary tag");
1011
1012 while ( adj && (adj != start) ) {
1013 pt = &mesh->tetra[adj];
1014
1015 /* identification of edge number in tetra adj */
1016 if ( !MMG3D_findEdge(mesh,pt,adj,na,nb,1,NULL,&i) ) return -1;
1017
1018 /* update edge ref and tag */
1019 if ( pt->xt ) {
1020 pxt = &mesh->xtetra[pt->xt];
1021
1022 /* if tag and edge are already consistent, no need to update the shell (do
1023 * not consider edges with tag 0 as they com from the creation of a xtetra
1024 * during a previous collapse and are not updated) */
1025 xtag = pxt->tag[i] | MG_BDY;
1026
1027 if ( pxt->tag[i] & MG_BDY ) {
1028 if ( xtag == tag && pxt->edg[i] == ref ) {
1029 return start;
1030 }
1031 }
1032
1033 pxt->edg[i] = ref;
1034 pxt->tag[i] |= tag;
1035 if( ((xtag & MG_REQ) && !(xtag & MG_NOSURF)) ||
1036 (( tag & MG_REQ) && !( tag & MG_NOSURF)))
1037 pxt->tag[i] &= ~MG_NOSURF;
1038 }
1039
1040 /* set new triangle for travel */
1041 adja = &mesh->adja[4*(adj-1)+1];
1042 if ( pt->v[ MMG5_ifar[i][0] ] == piv ) {
1043 adj = adja[ MMG5_ifar[i][0] ] / 4;
1044 piv = pt->v[ MMG5_ifar[i][1] ];
1045 }
1046 else {
1047 adj = adja[ MMG5_ifar[i][1] ] /4;
1048 piv = pt->v[ MMG5_ifar[i][0] ];
1049 }
1050 }
1051
1052 return adj;
1053}
1054
1066static inline
1067int MMG3D_update_shellEdgeTag(MMG5_pMesh mesh,MMG5_int start, int8_t ia,int16_t tag,MMG5_int ref) {
1068 MMG5_pTetra pt;
1069 MMG5_pxTetra pxt;
1070 MMG5_int piv,na,nb,adj,*adja;
1071 int16_t xtag;
1072
1073 pt = &mesh->tetra[start];
1074
1075 assert( start >= 1 && MG_EOK(pt) && "invalid tetra" );
1076 assert ( tag & MG_BDY && "Unexpected non boundary tag");
1077
1078 pxt = NULL;
1079 na = pt->v[MMG5_iare[ia][0]];
1080 nb = pt->v[MMG5_iare[ia][1]];
1081
1082 if ( pt->xt ) {
1083 pxt = &mesh->xtetra[pt->xt];
1084
1085 /* if tag and edge are already consistent, no need to update the shell (do
1086 * not consider edges with tag 0 as they com from the creation of a xtetra
1087 * during a previous collapse and are not updated) */
1088 xtag = pxt->tag[ia] | MG_BDY;
1089
1090 if ( pxt->tag[ia] & MG_BDY ) {
1091 if ( xtag == tag && pxt->edg[ia] == ref ) {
1092 return 1;
1093 }
1094 }
1095
1096 pxt->tag[ia] |= tag;
1097 if( ((xtag & MG_REQ) && !(xtag & MG_NOSURF)) ||
1098 (( tag & MG_REQ) && !( tag & MG_NOSURF)))
1099 pxt->tag[ia] &= ~MG_NOSURF;
1100 pxt->edg[ia] = ref;
1101 }
1102
1103 /* Travel in one direction */
1104 adja = &mesh->adja[4*(start-1)+1];
1105 adj = adja[MMG5_ifar[ia][0]] / 4;
1106 piv = pt->v[MMG5_ifar[ia][1]];
1107
1108 adj = MMG3D_update_shellEdgeTag_oneDir(mesh,start,na,nb,tag,ref,piv,adj);
1109
1110 /* If all shell has been travelled, stop, else, travel it the other sense */
1111 if ( adj == start ) return 1;
1112 else if ( adj < 0 ) return 0;
1113
1114 assert(!adj);
1115
1116 pt = &mesh->tetra[start];
1117 adja = &mesh->adja[4*(start-1)+1];
1118 adj = adja[MMG5_ifar[ia][1]] / 4;
1119 piv = pt->v[MMG5_ifar[ia][0]];
1120
1121 adj = MMG3D_update_shellEdgeTag_oneDir(mesh,start,na,nb,tag,ref,piv,adj);
1122
1123 if ( adj < 0 ) return 0;
1124
1125 return 1;
1126}
1127
1144MMG5_int MMG5_colver(MMG5_pMesh mesh,MMG5_pSol met,int64_t *list,int ilist,int8_t indq,int8_t typchk) {
1145 MMG5_pTetra pt,pt1;
1146 MMG5_pxTetra pxt,pxt1;
1147 MMG5_xTetra xt,xts;
1148 MMG5_int iel,jel,pel,qel,k,np,nq,*adja;
1149 uint8_t ip,iq,j,voy,voyp,voyq;
1150
1151 /* coledge[i] contains the local indices of edges that will be merged by the
1152 * collapse corresponding with the configuration i. The edge coledge[i][0] is
1153 * merged with the edge coledge[i][1] a,d the edge coledge[i][2] is merged
1154 * with edge coledge[i][3].
1155 * Config 0: merge of vertices 0 and 1
1156 * config 1: merge of vertices 0 and 2
1157 * config 2: merge of vertices 0 and 3
1158 * config 3: merge of vertices 1 and 2
1159 * config 4: merge of vertices 1 and 3
1160 * config 5: merge of vertices 2 and 3
1161 */
1162 const int8_t MMG5_coledge[6][4] = {
1163 {1,3,2,4}, {0,3,2,5}, {0,4,1,5},{0,1,4,5}, {0,2,3,5}, {4,3,2,1} };
1164
1165 iel = list[0] / 4;
1166 ip = list[0] % 4;
1167 pt = &mesh->tetra[iel];
1168 np = pt->v[ip];
1169 nq = pt->v[indq];
1170
1171 /* Mark elements of the shell of edge (pq) */
1172 for (k=0; k<ilist; k++) {
1173 iel = list[k] / 4;
1174 ip = list[k] % 4;
1175 pt = &mesh->tetra[iel];
1176
1177 for (j=0; j<3; j++) {
1178 iq = MMG5_idir[ip][j];
1179 if ( pt->v[iq] == nq ) {
1180
1181 list[k] *= -1;
1182 break;
1183 }
1184 }
1185 }
1186
1187 /* Avoid recreating existing elt */
1188 for (k=0; k<ilist; k++) {
1189
1190 if ( list[k] < 0 ) continue;
1191
1192 iel = list[k] / 4;
1193 ip = list[k] % 4;
1194
1195 adja = &mesh->adja[4*(iel-1)+1];
1196 jel = adja[ip];
1197 if ( !jel ) continue;
1198
1199 jel /= 4;
1200 voy = adja[ip] % 4;
1201 pt = &mesh->tetra[jel];
1202 if (pt->v[voy] == nq) {
1203 return 0;
1204 }
1205 }
1206
1212 for (k=0; k<ilist; k++) {
1213
1214 if ( list[k] > 0 ) continue;
1215
1216 iel = -list[k] / 4;
1217 ip = -list[k] % 4;
1218 pt = &mesh->tetra[iel];
1219
1220 for (j=0; j<3; j++) {
1221 iq = MMG5_idir[ip][j];
1222 if ( pt->v[iq] == nq ) {
1223 break;
1224 }
1225 }
1226
1227 /* The configuration is computed by setting the ip and iq bits to 1 */
1228 int8_t cfg = 0;
1229 MG_SET(cfg,ip);
1230 MG_SET(cfg,iq);
1231
1232 const int8_t *coled;
1233 switch(cfg) {
1234 case 3:
1235 /* collapse of vertices 0 and 1 */
1236 coled = MMG5_coledge[0];
1237 break;
1238 case 5:
1239 /* collapse of vertices 0 and 2 */
1240 coled = MMG5_coledge[1];
1241 break;
1242 case 9:
1243 /* collapse of vertices 0 and 3 */
1244 coled = MMG5_coledge[2];
1245 break;
1246 case 6:
1247 /* collapse of vertices 1 and 2 */
1248 coled = MMG5_coledge[3];
1249 break;
1250 case 10:
1251 /* collapse of vertices 1 and 3 */
1252 coled = MMG5_coledge[4];
1253 break;
1254 case 12:
1255 /* collapse of vertices 2 and 3 */
1256 coled = MMG5_coledge[5];
1257 break;
1258 default:
1259 assert ( 0 && "Unexpected collapse configuration");
1260 }
1261
1262 int l = 0;
1263 for ( l=0; l<3; l+=2 ) {
1264 /* when l=0 we update 2 edges that will be merged together, when l=2 we
1265 * update the two others. In practice we don't care which edge comes
1266 * from ip and which one from iq, it only matters that iped and iqed
1267 * will be merged at the end of the collapse. */
1268 int iped = coled[l+0];
1269 int iqed = coled[l+1];
1270
1271 int16_t tagip = 0;
1272 MMG5_int refip = 0;
1273 int16_t tagiq = 0;
1274 MMG5_int refiq = 0;
1275
1276 if ( !MMG3D_get_shellEdgeTag(mesh,iel,iped,&tagip,&refip) ) {
1277 fprintf(stderr,"\n ## Error: %s: 1. unable to get edge info.\n",__func__);
1278 return 0;
1279 }
1280 if ( !MMG3D_get_shellEdgeTag(mesh,iel,iqed,&tagiq,&refiq) ) {
1281 fprintf(stderr,"\n ## Error: %s: 2. unable to get edge info.\n",__func__);
1282 return 0;
1283 }
1284
1285 if ( (tagip != tagiq) || (refip != refiq) ) {
1286 /* Update edge infos */
1287 tagip |= tagiq;
1288 refip = (refip > refiq )? refip : refiq;
1289
1290 /* If the xtetra info are consistent when entering in colver taged and
1291 * refed contains suitable values and we can travel the shell of
1292 * iped and iqed to update all the needed info on edges */
1293 if ( !MMG3D_update_shellEdgeTag(mesh,iel,iped,tagip,refip) ) {
1294 fprintf(stderr,"\n ## Error: %s: 1. unable to update edge info.\n",__func__);
1295 return 0;
1296 }
1297 if ( !MMG3D_update_shellEdgeTag(mesh,iel,iqed,tagip,refip) ) {
1298 fprintf(stderr,"\n ## Error: %s: 1. unable to update edge info.\n",__func__);
1299 return 0;
1300 }
1301 }
1302#ifndef NDEBUG
1303 else {
1304 assert ( tagip == tagiq && refip == refiq );
1305 if ( mesh->info.ddebug && (tagip & MG_BDY) ) {
1306 MMG3D_chk_shellEdgeTag(mesh,iel,iped,tagip,refip);
1307 MMG3D_chk_shellEdgeTag(mesh,iel,iqed,tagip,refip);
1308 }
1309 }
1310#endif
1311 }
1312 }
1313
1314 /* deal with the shell of edge (pq) and the implied updates */
1315 for (k=0; k<ilist; k++) {
1316 if ( list[k] > 0 ) continue;
1317 iel = (-list[k]) / 4;
1318 ip = (-list[k]) % 4;
1319 pt = &mesh->tetra[iel];
1320
1321 iq = MMG5_inxt3[ip];
1322 for (j=0; j<3; j++) {
1323 if ( pt->v[iq] == nq ) break;
1324 iq = MMG5_inxt3[iq];
1325 }
1326 assert(j<3);
1327
1328
1329 adja = &mesh->adja[4*(iel-1)+1];
1330
1331 /* pel = neighbour of iel that belongs to ball of p \setminus shell, same for qel */
1332 pel = adja[iq] / 4;
1333 voyp = adja[iq] % 4;
1334 qel = adja[ip] / 4;
1335 voyq = adja[ip] % 4;
1336
1337 /* Update adjacency relations */
1338 if ( pel ) {
1339 adja = &mesh->adja[4*(pel-1)+1];
1340 adja[voyp] = 4*qel+voyq;
1341 }
1342 if ( qel ) {
1343 adja = &mesh->adja[4*(qel-1)+1];
1344 adja[voyq] = 4*pel+voyp;
1345 }
1346
1347 /* Update references for faces (one in pel) ;
1348 possibly, creation of a new field pxt for pel must be carried out */
1349 if ( pel ) {
1350 pt1 = &mesh->tetra[pel];
1351 if ( pt->xt ) {
1352 pxt = &mesh->xtetra[pt->xt];
1353 memcpy(&xts,pxt,sizeof(MMG5_xTetra));
1354 if ( pt1->xt > 0 ) {
1355 pxt1 = &mesh->xtetra[pt1->xt];
1356 pxt1->ref[voyp] = MG_MAX(pxt1->ref[voyp],pxt->ref[ip]);
1357 pxt1->ftag[voyp] = pxt1->ftag[voyp] | pxt->ftag[ip];
1358
1359 /* the face voyp of pxt1 is projected into the face ip of pxt so it
1360 * takes its orientation */
1361 if ( pxt->ftag[ip] & MG_BDY ) {
1362 if ( MG_GET(pxt->ori,ip) )
1363 MG_SET( pxt1->ori,voyp );
1364 else
1365 MG_CLR( pxt1->ori,voyp );
1366 }
1367 /* correct the orientation if the new adjacent has a different
1368 * reference */
1369 if( qel ) {
1370 if( pt1->ref > mesh->tetra[qel].ref )
1371 MG_SET( pxt1->ori,voyp );
1372 else if( pt1->ref < mesh->tetra[qel].ref )
1373 MG_CLR( pxt1->ori,voyp );
1374 }
1375
1376#ifndef NDEBUG
1377 else {
1378 /* Check that a non parallel external boundary face has always a good
1379 * orientation */
1380 if ( (pxt1->ftag[voyp] & MG_BDY) && !(pxt1->ftag[voyp] & MG_PARBDY) ) {
1381 assert ( MG_GET( pxt1->ori,voyp ) );
1382 }
1383 }
1384#endif
1385
1386 /* update tags for edges */
1387 MMG3D_update_edgeTag(pt,pxt,np,nq,ip,pt1,pxt1,voyp);
1388 }
1389 else {
1390 /* shell tet has a xtetra and pel exists but pel don't have a xtetra */
1391 pxt1 = &xt;
1392 memset(pxt1,0,sizeof(MMG5_xTetra));
1393 pxt1->ref[voyp] = pxt->ref[ip];
1394 pxt1->ftag[voyp] = pxt->ftag[ip];
1395 pxt1->ori = 15;
1396 if ( !MG_GET(pxt->ori,ip) ) MG_CLR(pxt1->ori,voyp);
1397
1398 /* update tags for edges */
1399 MMG3D_update_edgeTag(pt,pxt,np,nq,ip,pt1,pxt1,voyp);
1400
1401 /* Recover the already used place by pxt: now pel has a xtetra but the
1402 * edges coming from voyp (not directly involved in the collapse) will
1403 * have the tag 0 that may be non consistent with other tags in their
1404 * respective shells. */
1405 pt1->xt = pt->xt;
1406 memcpy(pxt,pxt1,sizeof(MMG5_xTetra));
1407 }
1408 }
1409 else {
1410 /* Shell tet don't have a xtetra: only the values of pel corresponding
1411 * to pt become 0 */
1412 if ( pt1->xt > 0 ) {
1413 pxt1 = &mesh->xtetra[pt1->xt];
1414 pxt1->ref[voyp] = 0;
1415 pxt1->ftag[voyp] = 0;
1416 MG_SET(pxt1->ori,voyp);
1417 }
1418 }
1419
1420 if ( qel ) {
1421 pt1 = &mesh->tetra[qel];
1422 if ( pt->xt ) {
1423 pxt = &xts;
1424 if ( pt1->xt > 0 ) {
1425 pxt1 = &mesh->xtetra[pt1->xt];
1426 pxt1->ref[voyq] = MG_MAX(pxt1->ref[voyq],pxt->ref[iq]);
1427 pxt1->ftag[voyq] = (pxt1->ftag[voyq] | pxt->ftag[iq]);
1428
1429 /* the face voyq of pxt1 is projected into the face iq of pxt so it
1430 * takes its orientation */
1431 if ( pxt->ftag[iq] & MG_BDY ) {
1432 if ( MG_GET(pxt->ori,iq) )
1433 MG_SET( pxt1->ori,voyq );
1434 else
1435 MG_CLR( pxt1->ori,voyq );
1436 }
1437 /* correct the orientation if the new adjacent has a different
1438 * reference */
1439 if( pel ) {
1440 if( pt1->ref > mesh->tetra[pel].ref )
1441 MG_SET( pxt1->ori,voyq );
1442 else if( pt1->ref < mesh->tetra[pel].ref )
1443 MG_CLR( pxt1->ori,voyq );
1444 }
1445
1446#ifndef NDEBUG
1447 else {
1448 /* Check that a non parallel external boundary face has always a good
1449 * orientation */
1450 if ( (pxt1->ftag[voyq] & MG_BDY) && !(pxt1->ftag[voyq] & MG_PARBDY) ) {
1451 assert ( MG_GET( pxt1->ori,voyq ) );
1452 }
1453 }
1454#endif
1455
1456 /* update tags for edges */
1457 MMG3D_update_edgeTag(pt,pxt,nq,np,iq,pt1,pxt1,voyq);
1458 }
1459 else {
1460 /* pel exists, shell tet has a xtetra but qel doesn't have boundary
1461 * tetra: create it */
1462 pxt1 = &xt;
1463 memset(pxt1,0,sizeof(MMG5_xTetra));
1464 pxt1->ref[voyq] = pxt->ref[iq];
1465 pxt1->ftag[voyq] = pxt->ftag[iq];
1466 pxt1->ori = 15;
1467 if ( !MG_GET(pxt->ori,iq) ) MG_CLR(pxt1->ori,voyq);
1468
1469 /* update tags for edges */
1470 MMG3D_update_edgeTag(pt,pxt,nq,np,iq,pt1,pxt1,voyq);
1471
1472 /* Create new field xt */
1473 mesh->xt++;
1474 if ( mesh->xt > mesh->xtmax ) {
1476 "larger xtetra table",
1477 mesh->xt--;return -1;);
1478 }
1479 /* Now qel has a xtetra but the edges coming from voyq (not directly
1480 * involved in the collapse) will have the tag 0 that may be non
1481 * consistent with other tags in their respective shells. */
1482 pt1->xt = mesh->xt;
1483 pxt = &mesh->xtetra[pt1->xt];
1484 memcpy(pxt,pxt1,sizeof(MMG5_xTetra));
1485 }
1486 }
1487 else {
1488 /* pel exist but shell tet doesn't have a boundary tetra: only the
1489 * values of qel corresponding to pt become 0 */
1490 if ( pt1->xt > 0 ) {
1491 pxt1 = &mesh->xtetra[pt1->xt];
1492 pxt1->ref[voyq] = 0;
1493 pxt1->ftag[voyq] = 0;
1494 MG_SET(pxt1->ori,voyq);
1495 }
1496 }
1497 }
1498 }
1499 else {
1500 /* pel==0: No adjacent through face iq */
1501 assert(pt->xt);
1502 pxt = &mesh->xtetra[pt->xt];
1503 if ( qel ) {
1504 pt1 = &mesh->tetra[qel];
1505 if ( pt1->xt > 0 ) {
1506 pxt1 = &mesh->xtetra[pt1->xt];
1507 pxt1->ref[voyq] = pxt->ref[iq];
1508 pxt1->ftag[voyq] = pxt->ftag[iq];
1509
1510 MG_SET(pxt1->ori,voyq);
1511
1512 /* update tags for edges */
1513 MMG3D_update_edgeTag(pt,pxt,nq,np,iq,pt1,pxt1,voyq);
1514 }
1515 else {
1516 pxt1 = &xt;
1517 memset(pxt1,0,sizeof(MMG5_xTetra));
1518 pxt1->ref[voyq] = pxt->ref[iq];
1519 pxt1->ftag[voyq] = pxt->ftag[iq];
1520 pxt1->ori = 15;
1521
1522 /* update tags for edges */
1523 MMG3D_update_edgeTag(pt,pxt,nq,np,iq,pt1,pxt1,voyq);
1524
1525 /* Recover the already used place by pxt */
1526 pt1->xt = pt->xt;
1527 memcpy(pxt,pxt1,sizeof(MMG5_xTetra));
1528 }
1529 }
1530 }
1531 if ( !MMG3D_delElt(mesh,iel) ) {
1532 return -1;
1533 }
1534 }
1535
1536 /* Update vertices coordinates for elements that do not belong to the shell of (pq) */
1537 for (k=0; k<ilist; k++) {
1538 if ( list[k] < 0 ) continue;
1539 iel = list[k] / 4;
1540 ip = list[k] % 4;
1541 pt = &mesh->tetra[iel];
1542 pt->v[ip] = nq;
1543 if ( typchk==1 && met->m && met->size > 1 )
1544 pt->qual=MMG5_caltet33_ani(mesh,met,pt);
1545 else
1546 pt->qual=MMG5_orcal(mesh,met,iel);
1547 pt->mark=mesh->mark;
1548 }
1549
1550 return np;
1551}
int ier
MMG5_pMesh * mesh
MMG5_Info info
int MMG3D_findEdge(MMG5_pMesh mesh, MMG5_pTetra pt, MMG5_int k, MMG5_int na, MMG5_int nb, int error, int8_t *mmgWarn, int8_t *ia)
Definition: boulep_3d.c:113
int MMG3D_chk_shellEdgeTag(MMG5_pMesh mesh, MMG5_int start, int8_t ia, int16_t tag, MMG5_int ref)
Definition: chkmsh_3d.c:142
static void MMG3D_update_edgeTag(MMG5_pTetra pt, MMG5_pxTetra pxt, MMG5_int np, MMG5_int nq, uint8_t ip, MMG5_pTetra pt1, MMG5_pxTetra pxt1, uint8_t voyp)
Definition: colver_3d.c:940
static int MMG3D_get_shellEdgeTag_oneDir(MMG5_pMesh mesh, MMG5_int start, MMG5_int na, MMG5_int nb, int16_t *tag, MMG5_int *ref, MMG5_int piv, MMG5_int adj, int8_t *filled)
Definition: colver_3d.c:377
static int MMG5_topchkcol_bdy(MMG5_pMesh mesh, MMG5_int k, int iface, int8_t iedg, MMG5_int *lists, int ilists)
Definition: colver_3d.c:251
static int MMG3D_get_shellEdgeTag(MMG5_pMesh mesh, MMG5_int start, int8_t ia, int16_t *tag, MMG5_int *ref)
Definition: colver_3d.c:432
MMG5_int MMG5_colver(MMG5_pMesh mesh, MMG5_pSol met, int64_t *list, int ilist, int8_t indq, int8_t typchk)
Definition: colver_3d.c:1144
int8_t ddb
Definition: mmg3d1_delone.c:42
int MMG5_chkcol_int(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int k, int8_t iface, int8_t iedg, int64_t *list, int ilist, int8_t typchk)
Definition: colver_3d.c:43
static MMG5_int MMG3D_update_shellEdgeTag_oneDir(MMG5_pMesh mesh, MMG5_int start, MMG5_int na, MMG5_int nb, int16_t tag, MMG5_int ref, MMG5_int piv, MMG5_int adj)
Definition: colver_3d.c:1002
int MMG5_chkcol_bdy(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int k, int8_t iface, int8_t iedg, int64_t *listv, int ilistv, MMG5_int *lists, int ilists, MMG5_int refmin, MMG5_int refplus, int8_t typchk, int isnm, int8_t isnmint)
Definition: colver_3d.c:515
static int MMG3D_update_shellEdgeTag(MMG5_pMesh mesh, MMG5_int start, int8_t ia, int16_t tag, MMG5_int ref)
Definition: colver_3d.c:1067
static MMG5_int MMG3D_unfold_shell(MMG5_pMesh mesh, MMG5_int start, MMG5_int end, MMG5_int na, MMG5_int nb, MMG5_int piv, MMG5_int *iel, int8_t *iopp)
Definition: colver_3d.c:193
static double MMG5_orcal(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int iel)
static const int8_t MMG5_idirinv[4][4]
int MMG3D_normalAdjaTri(MMG5_pMesh, MMG5_int, int8_t, int, double n[3])
Definition: split_3d.c:466
int MMG5_startedgsurfball(MMG5_pMesh mesh, MMG5_int nump, MMG5_int numq, MMG5_int *list, int ilist)
Definition: tools_3d.c:91
int MMG3D_chkmanicoll(MMG5_pMesh, MMG5_int, int, int, MMG5_int, MMG5_int, MMG5_int, MMG5_int, int8_t, int8_t)
Definition: mmg3d2.c:1672
int MMG3D_delElt(MMG5_pMesh mesh, MMG5_int iel)
Definition: zaldy_3d.c:122
int8_t MMG5_chkedg(MMG5_pMesh mesh, MMG5_Tria *pt, int8_t ori, double, double, int)
Definition: mmg3d1.c:363
static const int8_t MMG5_iarf[4][3]
iarf[i]: edges of face opposite to vertex i
double MMG5_caltet33_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTetra pt)
Definition: quality_3d.c:109
static const uint8_t MMG5_iare[6][2]
vertices of extremities of the edges of the tetra
static const uint8_t MMG5_ifar[6][2]
ifar[i][]: faces sharing the ith edge of the tetra
void MMG5_tet2tri(MMG5_pMesh mesh, MMG5_int k, int8_t ie, MMG5_Tria *ptt)
Definition: mmg3d1.c:102
static const uint8_t MMG5_inxt3[7]
next vertex of tetra: {1,2,3,0,1,2,3}
int MMG5_norface(MMG5_pMesh mesh, MMG5_int k, int iface, double v[3])
Definition: tools_3d.c:52
static const uint8_t MMG5_idir[4][3]
idir[i]: vertices of face opposite to vertex i
#define MG_PLUS
Definition: libmmgtypes.h:71
@ MMG5_Tetrahedron
Definition: libmmgtypes.h:227
@ MMG5_Triangle
Definition: libmmgtypes.h:226
#define MG_MINUS
Definition: libmmgtypes.h:76
#define MG_Tria
#define MG_REQ
#define MG_GEO
#define MMG5_EPSOK
#define MG_EOK(pt)
#define MG_CLR(flag, bit)
#define MG_PARBDY
#define MG_MIN(a, b)
#define MG_MAX(a, b)
#define MMG5_GAP
static const uint8_t MMG5_iprv2[3]
#define MMG5_TAB_RECALLOC(mesh, ptr, initSize, wantedGap, type, message, law)
#define MMG5_NULKAL
#define MG_GEO_OR_NOM(tag)
#define MG_Tetra
#define MG_GET(flag, bit)
static const uint8_t MMG5_inxt2[6]
#define MG_BDY
#define MMG5_ANGEDG
#define MG_NOSURF
#define MG_SIN_OR_NOM(tag)
int MMG5_norpts(MMG5_pMesh, MMG5_int, MMG5_int, MMG5_int, double *)
Definition: tools.c:183
#define MG_SET(flag, bit)
int8_t iso
Definition: libmmgtypes.h:534
int8_t parTyp
Definition: libmmgtypes.h:541
int8_t ddebug
Definition: libmmgtypes.h:532
double hmax
Definition: libmmgtypes.h:518
MMG5_pPar par
Definition: libmmgtypes.h:517
double dhd
Definition: libmmgtypes.h:518
double hausd
Definition: libmmgtypes.h:518
int8_t fem
Definition: libmmgtypes.h:539
MMG mesh structure.
Definition: libmmgtypes.h:605
MMG5_Info info
Definition: libmmgtypes.h:651
MMG5_int xt
Definition: libmmgtypes.h:620
MMG5_pPoint point
Definition: libmmgtypes.h:641
MMG5_int mark
Definition: libmmgtypes.h:618
MMG5_int * adja
Definition: libmmgtypes.h:624
MMG5_int xtmax
Definition: libmmgtypes.h:612
MMG5_pTetra tetra
Definition: libmmgtypes.h:643
MMG5_pxTetra xtetra
Definition: libmmgtypes.h:644
double hmax
Definition: libmmgtypes.h:259
double hausd
Definition: libmmgtypes.h:260
MMG5_int ref
Definition: libmmgtypes.h:261
int8_t elt
Definition: libmmgtypes.h:262
Structure to store points of a MMG mesh.
Definition: libmmgtypes.h:270
int16_t tag
Definition: libmmgtypes.h:284
MMG5_int xp
Definition: libmmgtypes.h:279
double * m
Definition: libmmgtypes.h:671
MMG5_int v[4]
Definition: libmmgtypes.h:403
MMG5_int xt
Definition: libmmgtypes.h:407
MMG5_int ref
Definition: libmmgtypes.h:404
MMG5_int mark
Definition: libmmgtypes.h:406
double qual
Definition: libmmgtypes.h:402
int16_t tag[3]
Definition: libmmgtypes.h:342
MMG5_int ref
Definition: libmmgtypes.h:335
MMG5_int v[3]
Definition: libmmgtypes.h:334
Structure to store the surface tetrahedra of a MMG mesh.
Definition: libmmgtypes.h:418
int16_t ftag[4]
Definition: libmmgtypes.h:423
int16_t tag[6]
Definition: libmmgtypes.h:425
int8_t ori
Definition: libmmgtypes.h:427
MMG5_int edg[6]
Definition: libmmgtypes.h:421
MMG5_int ref[4]
Definition: libmmgtypes.h:419