Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
movpt_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
36#include "libmmg3d.h"
39#include "mmgexterns_private.h"
40
59 int64_t *list,int ilist,int improve) {
60 MMG5_pTetra pt,pt0;
61 MMG5_pPoint p0,p1,p2,p3,ppt0;
62 double vol,totvol;
63 double calold,calnew,callist[MMG3D_LMAX+2];
64 double len1,len2;
65 int iloc;
66 int k,i0;
67 MMG5_int iel;
68
69 pt0 = &mesh->tetra[0];
70 ppt0 = &mesh->point[0];
71 memset(ppt0,0,sizeof(MMG5_Point));
72
73 /* Coordinates of optimal point */
74 calold = DBL_MAX;
75 totvol = 0.0;
76 for (k=0; k<ilist; k++) {
77 iel = list[k] / 4;
78 pt = &mesh->tetra[iel];
79 p0 = &mesh->point[pt->v[0]];
80 p1 = &mesh->point[pt->v[1]];
81 p2 = &mesh->point[pt->v[2]];
82 p3 = &mesh->point[pt->v[3]];
83 vol= MMG5_det4pt(p0->c,p1->c,p2->c,p3->c);
84 totvol += vol;
85 /* barycenter */
86 ppt0->c[0] += 0.25 * vol*(p0->c[0] + p1->c[0] + p2->c[0] + p3->c[0]);
87 ppt0->c[1] += 0.25 * vol*(p0->c[1] + p1->c[1] + p2->c[1] + p3->c[1]);
88 ppt0->c[2] += 0.25 * vol*(p0->c[2] + p1->c[2] + p2->c[2] + p3->c[2]);
89 calold = MG_MIN(calold, pt->qual);
90 }
91 if (totvol < MMG5_EPSD2) {
92 return 0;
93 }
94
95 totvol = 1.0 / totvol;
96 ppt0->c[0] *= totvol;
97 ppt0->c[1] *= totvol;
98 ppt0->c[2] *= totvol;
99
100 /* Check new position validity */
101 calnew = DBL_MAX;
102 i0 = -1;
103 assert ( ilist>0 );
104 for (k=0; k<ilist; k++) {
105 iel = list[k] / 4;
106 i0 = list[k] % 4;
107 pt = &mesh->tetra[iel];
108 memcpy(pt0,pt,sizeof(MMG5_Tetra));
109 pt0->v[i0] = 0;
110 callist[k] = MMG5_orcal(mesh,met,0);
111 if (callist[k] < MMG5_NULKAL) {
112 return 0;
113 }
114 calnew = MG_MIN(calnew,callist[k]);
115
116 if ( improve==2 ) {
117 for ( iloc = 0; iloc < 3; ++iloc ) {
118 len1 = MMG5_lenedg_iso(mesh,met,MMG5_arpt[i0][iloc],pt);
119 len2 = MMG5_lenedg_iso(mesh,met,MMG5_arpt[i0][iloc],pt0);
120 if ( (len1 < MMG3D_LOPTL && len2 >= MMG3D_LOPTL) ||
121 (len1 > MMG3D_LOPTL && len2 >len1 ) ) {
122 return 0;
123 }
124
125 if ( (len1 > MMG3D_LOPTS && len2 <= MMG3D_LOPTS) ||
126 (len1 < MMG3D_LOPTS && len2 <len1 ) ) {
127 return 0;
128 }
129 }
130 }
131
132 }
133 if (calold < MMG5_EPSOK && calnew <= calold) {
134 return 0;
135 }
136 else if (calnew < MMG5_EPSOK) {
137 return 0;
138 }
139 else if ( improve && calnew < 1.02 * calold ) {
140 return 0;
141 }
142 else if ( calnew < 0.3 * calold ) {
143 return 0;
144 }
145
146 /* update position */
147 assert ( i0 >=0 );
148 if ( PROctree )
149 MMG3D_movePROctree(mesh, PROctree, pt->v[i0], ppt0->c, p0->c);
150
151 p0 = &mesh->point[pt->v[i0]];
152 p0->c[0] = ppt0->c[0];
153 p0->c[1] = ppt0->c[1];
154 p0->c[2] = ppt0->c[2];
155 for (k=0; k<ilist; k++) {
156 (&mesh->tetra[list[k]/4])->qual=callist[k];
157 (&mesh->tetra[list[k]/4])->mark=mesh->mark;
158 }
159
160 return 1;
161}
162
186 MMG5_int *list,int ilist,int improve) {
187 MMG5_pTetra pt,pt0;
188 MMG5_pPoint p0,p1,p2,p3,ppt0;
189 double vol,totvol;
190 double calold,calnew,callist[MMG3D_LMAX+2];
191 double x21,y21,z21,x31,y31,z31,nx,ny,nz,bary[3],dd,len;
192 double u10[3],u20[3],u30[3],oldc[3],coe;
193 int k,ifac,iter,maxtou;
194 MMG5_int iel;
195
196 pt0 = &mesh->tetra[0];
197 ppt0 = &mesh->point[0];
198 memset(ppt0,0,sizeof(MMG5_Point));
199
200 /* Coordinates of optimal point */
201 calold = DBL_MAX;
202 totvol = 0.0;
203 assert ( ilist );
204 for (k=0; k<ilist; k++) {
205 iel = list[k] / 4;
206 ifac = list[k] % 4;
207 pt = &mesh->tetra[iel];
208 assert ( MG_EOK(pt) );
209
210 p0 = &mesh->point[pt->v[ifac]];
211 memcpy(oldc,p0->c,3*sizeof(double));
212
213 p1 = &mesh->point[pt->v[MMG5_idir[ifac][0]]];
214 p2 = &mesh->point[pt->v[MMG5_idir[ifac][1]]];
215 p3 = &mesh->point[pt->v[MMG5_idir[ifac][2]]];
216
217 vol= MMG5_det4pt(p0->c,p1->c,p2->c,p3->c);
218 totvol += vol;
219
220 /* local optimal point : length 1 over the normal of the face */
221 x21 = p2->c[0]-p1->c[0];
222 y21 = p2->c[1]-p1->c[1];
223 z21 = p2->c[2]-p1->c[2];
224
225 x31 = p3->c[0]-p1->c[0];
226 y31 = p3->c[1]-p1->c[1];
227 z31 = p3->c[2]-p1->c[2];
228
229 nx = (y31*z21 - z31*y21);
230 ny = (z31*x21 - x31*z21);
231 nz = (x31*y21 - y31*x21);
232
233 dd = sqrt(nx*nx+ny*ny+nz*nz);
234 assert ( dd > 0. && "degenerated element");
235 dd = 1./dd;
236 nx *= dd;
237 ny *= dd;
238 nz *= dd;
239
240 u10[0] = p1->c[0]-p0->c[0];
241 u10[1] = p1->c[1]-p0->c[1];
242 u10[2] = p1->c[2]-p0->c[2];
243
244 u20[0] = p2->c[0]-p0->c[0];
245 u20[1] = p2->c[1]-p0->c[1];
246 u20[2] = p2->c[2]-p0->c[2];
247
248 u30[0] = p3->c[0]-p0->c[0];
249 u30[1] = p3->c[1]-p0->c[1];
250 u30[2] = p3->c[2]-p0->c[2];
251
252 assert ( met->m[pt->v[MMG5_idir[ifac][0]]] > 0. && "null metric at vertex 0" );
253 assert ( met->m[pt->v[MMG5_idir[ifac][1]]] > 0. && "null metric at vertex 1" );
254 assert ( met->m[pt->v[MMG5_idir[ifac][2]]] > 0. && "null metric at vertex 2" );
255
256 len = sqrt(u10[0]*u10[0]+u10[1]*u10[1]+u10[2]*u10[2])/met->m[pt->v[MMG5_idir[ifac][0]]]
257 + sqrt(u20[0]*u20[0]+u20[1]*u20[1]+u20[2]*u20[2])/met->m[pt->v[MMG5_idir[ifac][1]]]
258 + sqrt(u30[0]*u30[0]+u30[1]*u30[1]+u30[2]*u30[2])/met->m[pt->v[MMG5_idir[ifac][2]]];
259
260 len /= 3.;
261
262 assert ( len>0. && "degenerated element" );
263
264 len = 1./len;
265
266 /* face barycenter */
267 bary[0] = (p1->c[0]+p2->c[0]+p3->c[0])/3.;
268 bary[1] = (p1->c[1]+p2->c[1]+p3->c[1])/3.;
269 bary[2] = (p1->c[2]+p2->c[2]+p3->c[2])/3.;
270
271 /* Optimal point: barycenter of local optimal points */
272 vol = 1;
273 ppt0->c[0] += vol* ( bary[0] + nx * len );
274 ppt0->c[1] += vol* ( bary[1] + ny * len );
275 ppt0->c[2] += vol* ( bary[2] + nz * len );
276 calold = MG_MIN(calold, pt->qual);
277 }
278 if (totvol < MMG5_EPSD2) {
279 return 0;
280 }
281
282 ppt0->c[0] *= 1./(double) ilist;
283 ppt0->c[1] *= 1./(double) ilist;
284 ppt0->c[2] *= 1./(double) ilist;
285
286 coe = 0.9;
287 maxtou = 10;
288 iter = 0;
289 do {
290 p0->c[0] = (1. - coe) *oldc[0] + coe * ppt0->c[0];
291 p0->c[1] = (1. - coe) *oldc[1] + coe * ppt0->c[1];
292 p0->c[2] = (1. - coe) *oldc[2] + coe * ppt0->c[2];
293
294 /* Check new position validity */
295 calnew = DBL_MAX;
296 for (k=0; k<ilist; k++) {
297 iel = list[k] / 4;
298 pt = &mesh->tetra[iel];
299 memcpy(pt0,pt,sizeof(MMG5_Tetra));
300 callist[k] = MMG5_caltet(mesh,met,pt0); // MMG5_orcal(mesh,met,0)
301 if (calold < MMG5_EPSOK && callist[k] <= calold) {
302 break;
303 } else if ((callist[k] < MMG5_EPSOK)) {
304 break;
305 } else if( callist[k] < 0.1) {
306 if (callist[k] < 1.01*calold) break;
307 } else {
308 if (callist[k] < 1.02*calold) break;
309 }
310 calnew = MG_MIN(calnew,callist[k]);
311 }
312 if(k<ilist) {
313 coe *= 0.5;
314 continue;
315 }
316 break;
317 } while(++iter <= maxtou );
318
319 if ( iter > maxtou ) {
320 memcpy(p0->c,oldc,3*sizeof(double));
321 return 0;
322 }
323
324 /* update position */
325 if ( PROctree )
326 MMG3D_movePROctree(mesh, PROctree, mesh->tetra[list[0]/4].v[list[0]%4], ppt0->c, p0->c);
327
328 for (k=0; k<ilist; k++) {
329 (&mesh->tetra[list[k]/4])->qual=callist[k];
330 (&mesh->tetra[list[k]/4])->mark=mesh->mark;
331 }
332
333 return 1;
334}
335
350int MMG3D_rotate_surfacicBall(MMG5_pMesh mesh,MMG5_int *lists,int ilists,MMG5_int ip0,
351 double r[3][3],double *lispoi) {
352 MMG5_pTetra pt;
353 MMG5_pPoint p0,p1;
354 double ux,uy,uz,det2d;
355 MMG5_int k,na,nb,ntempa,ntempb;
356 int l;
357 uint8_t iface,i;
358
359 k = lists[0] / 4;
360 iface = lists[0] % 4;
361 pt = &mesh->tetra[k];
362 p0 = &mesh->point[ip0];
363
364 na = nb = 0;
365 for (i=0; i<3; i++) {
366 if ( pt->v[MMG5_idir[iface][i]] != ip0 ) {
367 if ( !na )
368 na = pt->v[MMG5_idir[iface][i]];
369 else
370 nb = pt->v[MMG5_idir[iface][i]];
371 }
372 }
373
374 for (l=1; l<ilists; l++) {
375 k = lists[l] / 4;
376 iface = lists[l] % 4;
377 pt = &mesh->tetra[k];
378 ntempa = ntempb = 0;
379 for (i=0; i<3; i++) {
380 if ( pt->v[MMG5_idir[iface][i]] != ip0 ) {
381 if ( !ntempa )
382 ntempa = pt->v[MMG5_idir[iface][i]];
383 else
384 ntempb = pt->v[MMG5_idir[iface][i]];
385 }
386 }
387 if ( ntempa == na )
388 p1 = &mesh->point[na];
389 else if ( ntempa == nb )
390 p1 = &mesh->point[nb];
391 else if ( ntempb == na )
392 p1 = &mesh->point[na];
393 else {
394 assert(ntempb == nb);
395 p1 = &mesh->point[nb];
396 }
397 ux = p1->c[0] - p0->c[0];
398 uy = p1->c[1] - p0->c[1];
399 uz = p1->c[2] - p0->c[2];
400
401 lispoi[3*l+1] = r[0][0]*ux + r[0][1]*uy + r[0][2]*uz;
402 lispoi[3*l+2] = r[1][0]*ux + r[1][1]*uy + r[1][2]*uz;
403 lispoi[3*l+3] = r[2][0]*ux + r[2][1]*uy + r[2][2]*uz;
404
405 na = ntempa;
406 nb = ntempb;
407 }
408
409 /* Finish with point 0 */;
410 k = lists[0] / 4;
411 iface = lists[0] % 4;
412 pt = &mesh->tetra[k];
413 ntempa = ntempb = 0;
414 for (i=0; i<3; i++) {
415 if ( pt->v[MMG5_idir[iface][i]] != ip0 ) {
416 if ( !ntempa )
417 ntempa = pt->v[MMG5_idir[iface][i]];
418 else
419 ntempb = pt->v[MMG5_idir[iface][i]];
420 }
421 }
422 if ( ntempa == na )
423 p1 = &mesh->point[na];
424 else if ( ntempa == nb )
425 p1 = &mesh->point[nb];
426 else if ( ntempb == na )
427 p1 = &mesh->point[na];
428 else {
429 assert(ntempb == nb);
430 p1 = &mesh->point[nb];
431 }
432
433 ux = p1->c[0] - p0->c[0];
434 uy = p1->c[1] - p0->c[1];
435 uz = p1->c[2] - p0->c[2];
436
437 lispoi[1] = r[0][0]*ux + r[0][1]*uy + r[0][2]*uz;
438 lispoi[2] = r[1][0]*ux + r[1][1]*uy + r[1][2]*uz;
439 lispoi[3] = r[2][0]*ux + r[2][1]*uy + r[2][2]*uz;
440
441 /* list goes modulo ilist */
442 lispoi[3*ilists+1] = lispoi[1];
443 lispoi[3*ilists+2] = lispoi[2];
444 lispoi[3*ilists+3] = lispoi[3];
445
449 /* Check all projections over tangent plane. */
450 for (k=0; k<ilists-1; k++) {
451 det2d = lispoi[3*k+1]*lispoi[3*(k+1)+2] - lispoi[3*k+2]*lispoi[3*(k+1)+1];
452 if ( det2d < 0.0 ) {
453 return 0;
454 }
455 }
456 det2d = lispoi[3*(ilists-1)+1]*lispoi[3*0+2] - lispoi[3*(ilists-1)+2]*lispoi[3*0+1];
457 if ( det2d < 0.0 ) {
458 return 0;
459 }
460
461 return 1;
462}
463
464
482int MMG3D_movbdyregpt_geom(MMG5_pMesh mesh,MMG5_int *lists,const MMG5_int kel,
483 const MMG5_int ip0,double n[3],double lambda[3],double o[3],
484 double no[3]) {
485 MMG5_pTetra pt;
486 MMG5_pxTetra pxt;
487 MMG5_pPoint p1,p2,ppt0,p0;
488 MMG5_Tria tt;
489 MMG5_pxPoint pxp;
490 MMG5_Bezier b;
491 double uv[2],to[3],detloc;
492 int iel,na,nb,ntempb,ntempc,nxp;
493 uint8_t iface,i;
494 static int8_t mmgErr0=0,mmgErr1=0;
495
496 iel = lists[kel] / 4;
497 iface = lists[kel] % 4;
498 pt = &mesh->tetra[iel];
499 pxt = &mesh->xtetra[pt->xt];
500 p0 = &mesh->point[ip0];
501
502 assert( 0<=iface && iface<4 && "unexpected local face idx");
503 MMG5_tet2tri(mesh,iel,iface,&tt);
504
505 if(!MMG5_bezierCP(mesh,&tt,&b,MG_GET(pxt->ori,iface))){
506 if( !mmgErr0 ) {
507 mmgErr0 = 1;
508 fprintf(stderr,"\n ## Error: %s: function MMG5_bezierCP return 0.\n",
509 __func__);
510 }
511 return -1;
512 }
513
514 /* Now, for Bezier interpolation, one should identify which of i,i1,i2 is
515 0,1,2 recall uv[0] = barycentric coord associated to pt->v[1], uv[1]
516 associated to pt->v[2] and 1-uv[0]-uv[1] is associated to pt->v[0]. For
517 this, use the fact that kel, kel + 1 is positively oriented with respect to
518 n */
519 na = nb = 0;
520 for( i=0 ; i<4 ; i++ ){
521 if ( (pt->v[i] != ip0) && (pt->v[i] != pt->v[iface]) ) {
522 if ( !na )
523 na = pt->v[i];
524 else
525 nb = pt->v[i];
526 }
527 }
528 p1 = &mesh->point[na];
529 p2 = &mesh->point[nb];
530 detloc = MMG5_det3pt1vec(p0->c,p1->c,p2->c,n);
531
532 /* ntempa=point to which is associated 1-uv[0]-uv[1], ntempb=uv[0], ntempc=uv[1] */
533 ntempb = pt->v[MMG5_idir[iface][1]];
534 ntempc = pt->v[MMG5_idir[iface][2]];
535
536 /* na = lispoi[kel] -> lambda[1], nb = lispoi[kel+1] -> lambda[2] */
537 if ( detloc > 0.0 ) {
538 if ( ntempb == na )
539 uv[0] = lambda[1];
540 else if (ntempb == nb )
541 uv[0] = lambda[2];
542 else {
543 assert(ntempb == ip0);
544 uv[0] = lambda[0];
545 }
546 if ( ntempc == na )
547 uv[1] = lambda[1];
548 else if (ntempc == nb )
549 uv[1] = lambda[2];
550 else {
551 assert(ntempc == ip0);
552 uv[1] = lambda[0];
553 }
554 }
555 /* nb = lispoi[kel] -> lambda[1], na = lispoi[kel+1] -> lambda[2] */
556 else {
557 if ( ntempb == na )
558 uv[0] = lambda[2];
559 else if (ntempb == nb )
560 uv[0] = lambda[1];
561 else {
562 assert(ntempb == ip0);
563 uv[0] = lambda[0];
564 }
565 if ( ntempc == na )
566 uv[1] = lambda[2];
567 else if (ntempc == nb )
568 uv[1] = lambda[1];
569 else {
570 assert(ntempc == ip0);
571 uv[1] = lambda[0];
572 }
573 }
574 if(!MMG3D_bezierInt(&b,uv,o,no,to)){
575 if( !mmgErr1 ) {
576 mmgErr1 = 1;
577 fprintf(stderr," ## Error: %s: function MMG3D_bezierInt return 0.\n",
578 __func__);
579 }
580 return -1;
581 }
582
583 /* Test : make sure that geometric approximation has not been degraded too much */
584 ppt0 = &mesh->point[0];
585 ppt0->c[0] = o[0];
586 ppt0->c[1] = o[1];
587 ppt0->c[2] = o[2];
588
589 ppt0->tag = p0->tag;
590 ppt0->ref = p0->ref;
591
592
593 nxp = mesh->xp + 1;
594 if ( nxp > mesh->xpmax ) {
596 "larger xpoint table",
597 return 0);
598 }
599 ppt0->xp = nxp;
600 pxp = &mesh->xpoint[nxp];
601 memcpy(pxp,&(mesh->xpoint[p0->xp]),sizeof(MMG5_xPoint));
602 pxp->n1[0] = no[0];
603 pxp->n1[1] = no[1];
604 pxp->n1[2] = no[2];
605
606 return nxp;
607}
608
626 int ilistv,MMG5_int *lists,int ilists,
627 int improveSurf,int improveVol) {
628 MMG5_pTetra pt,pt0;
629 MMG5_pPoint p0;
630 MMG5_Tria tt;
631 MMG5_pxPoint pxp;
632 double *n,r[3][3],lispoi[3*MMG3D_LMAX+1],ux,uy,det2d;
633 double detloc,oppt[2],step,lambda[3];
634 double ll,m[2],o[3],no[3];
635 double calold,calnew,caltmp,callist[MMG3D_LMAX+2];
636 int l,nut,nxp;
637 MMG5_int kel,k,ip0;
638 uint8_t i0,iface,i;
639
640 step = 0.1;
641 nut = 0;
642 oppt[0] = 0.0;
643 oppt[1] = 0.0;
644 if ( ilists < 2 ) return 0;
645
646 k = listv[0] / 4;
647 i0 = listv[0] % 4;
648 pt = &mesh->tetra[k];
649 ip0 = pt->v[i0];
650 p0 = &mesh->point[ip0];
651 assert( p0->xp && !MG_EDG(p0->tag) );
652
653 n = &(mesh->xpoint[p0->xp].n1[0]);
654
656 if ( !MMG5_rotmatrix(n,r) ) return 0;
657
660 if ( !MMG3D_rotate_surfacicBall(mesh,lists,ilists,ip0,r,lispoi) ) {
661 return 0;
662 }
663
666 for (k=0; k<ilists; k++) {
667 m[0] = 0.5*(lispoi[3*(k+1)+1] + lispoi[3*k+1]);
668 m[1] = 0.5*(lispoi[3*(k+1)+2] + lispoi[3*k+2]);
669 ux = lispoi[3*(k+1)+1] - lispoi[3*k+1];
670 uy = lispoi[3*(k+1)+2] - lispoi[3*k+2];
671 ll = ux*ux + uy*uy;
672 if ( ll < MMG5_EPSD2 ) continue;
673 nut++;
674 oppt[0] += (m[0]-MMG5_SQR32*uy);
675 oppt[1] += (m[1]+MMG5_SQR32*ux);
676 }
677 oppt[0] *= (1.0 / nut);
678 oppt[1] *= (1.0 / nut);
679
681 det2d = lispoi[1]*oppt[1] - lispoi[2]*oppt[0];
682 kel = 0;
683 if ( det2d >= 0.0 ) {
684 for (k=0; k<ilists; k++) {
685 detloc = oppt[0]*lispoi[3*(k+1)+2] - oppt[1]*lispoi[3*(k+1)+1];
686 if ( detloc >= 0.0 ) {
687 kel = k;
688 break;
689 }
690 }
691 if ( k == ilists ) return 0;
692 }
693 else {
694 for (k=ilists-1; k>=0; k--) {
695 detloc = lispoi[3*k+1]*oppt[1] - lispoi[3*k+2]*oppt[0];
696 if ( detloc >= 0.0 ) {
697 kel = k;
698 break;
699 }
700 }
701 if ( k == -1 ) return 0;
702 }
703
704 /* Sizing of time step : make sure point does not go out corresponding triangle. */
705 det2d = -oppt[1]*(lispoi[3*(kel+1)+1] - lispoi[3*(kel)+1]) + \
706 oppt[0]*(lispoi[3*(kel+1)+2] - lispoi[3*(kel)+2]);
707 if ( fabs(det2d) < MMG5_EPSD2 ) return 0;
708
709 det2d = 1.0 / det2d;
710 step *= det2d;
711
712 det2d = lispoi[3*(kel)+1]*(lispoi[3*(kel+1)+2] - lispoi[3*(kel)+2]) - \
713 lispoi[3*(kel)+2 ]*(lispoi[3*(kel+1)+1] - lispoi[3*(kel)+1]);
714 step *= det2d;
715 step = fabs(step);
716 oppt[0] *= step;
717 oppt[1] *= step;
718
719 /* Computation of the barycentric coordinates of the new point in the corresponding triangle. */
720 det2d = lispoi[3*kel+1]*lispoi[3*(kel+1)+2] - lispoi[3*kel+2]*lispoi[3*(kel+1)+1];
721 if ( det2d < MMG5_EPSD2 ) return 0;
722 det2d = 1.0 / det2d;
723 lambda[1] = lispoi[3*(kel+1)+2]*oppt[0] - lispoi[3*(kel+1)+1]*oppt[1];
724 lambda[2] = -lispoi[3*(kel)+2]*oppt[0] + lispoi[3*(kel)+1]*oppt[1];
725 lambda[1]*= (det2d);
726 lambda[2]*= (det2d);
727 lambda[0] = 1.0 - lambda[1] - lambda[2];
728
730 nxp = MMG3D_movbdyregpt_geom(mesh,lists,kel,ip0,n,lambda,o,no);
731 if ( nxp < 0 ) {
732 return -1;
733 }
734 else if ( !nxp ) {
735 return 0;
736 }
737 pxp = &mesh->xpoint[nxp];
738
739 /* For each surfacic triangle, build a virtual displaced triangle for check purposes */
740 calold = calnew = DBL_MAX;
741 for (l=0; l<ilists; l++) {
742 k = lists[l] / 4;
743 iface = lists[l] % 4;
744
745 assert( 0<=iface && iface<4 && "unexpected local face idx");
746 MMG5_tet2tri(mesh,k,iface,&tt);
747 calold = MG_MIN(calold,MMG5_caltri(mesh,met,&tt));
748
749 for( i=0 ; i<3 ; i++ )
750 if ( tt.v[i] == ip0 ) break;
751 assert(i<3);
752 if ( i==3 ) return 0;
753
754 tt.v[i] = 0;
755
756 caltmp = MMG5_caltri(mesh,met,&tt);
757 if ( caltmp < MMG5_EPSD2 ) {
758 /* We don't check the input triangle qualities, thus we may have a very
759 * bad triangle in our mesh */
760 return 0;
761 }
762 calnew = MG_MIN(calnew,caltmp);
763 }
764 if ( calold < MMG5_EPSOK && calnew <= calold ) return 0;
765 else if (calnew < MMG5_EPSOK) return 0;
766 else if (improveSurf && calnew < 1.02*calold) return 0;
767 else if ( calnew < 0.3*calold ) return 0;
768 memset(pxp,0,sizeof(MMG5_xPoint));
769
770 /* Test : check whether all volumes remain positive with new position of the point */
771 calold = calnew = DBL_MAX;
772
773 for (l=0; l<ilistv; l++) {
774 k = listv[l] / 4;
775 i0 = listv[l] % 4;
776 pt = &mesh->tetra[k];
777 pt0 = &mesh->tetra[0];
778 memcpy(pt0,pt,sizeof(MMG5_Tetra));
779 pt0->v[i0] = 0;
780 calold = MG_MIN(calold, pt->qual);
781 callist[l]=MMG5_orcal(mesh,met,0);
782
783 if (callist[l] < MMG5_NULKAL) {
784 return 0;
785 }
786 calnew = MG_MIN(calnew,callist[l]);
787
788 }
789
790
791 if (calold < MMG5_EPSOK && calnew <= calold) {
792 return 0;
793 }
794 else if (calnew < MMG5_EPSOK) {
795 return 0;
796 }
797 else if (improveVol && calnew < calold) {
798 return 0;
799 }
800 else if (calnew < 0.3*calold) {
801 return 0;
802 }
803
804 /* When all tests have been carried out, update coordinates and normals */
805 if ( PROctree )
806 MMG3D_movePROctree(mesh, PROctree, ip0, o, p0->c);
807
808 p0->c[0] = o[0];
809 p0->c[1] = o[1];
810 p0->c[2] = o[2];
811
812 n[0] = no[0];
813 n[1] = no[1];
814 n[2] = no[2];
815
816 for(l=0; l<ilistv; l++){
817 (&mesh->tetra[listv[l]/4])->qual= callist[l];
818 (&mesh->tetra[listv[l]/4])->mark=mesh->mark;
819 }
820 return 1;
821}
822
844static inline
845int MMG3D_curveEndingPts_chkEdg(MMG5_pMesh mesh,MMG5_int *lists,int l,MMG5_int ip0,
846 MMG5_int *ipa,MMG5_int *ipb,const int16_t edgTag,MMG5_int *ip) {
847
848 MMG5_pTetra pt;
849 MMG5_int iel,iptmpa,iptmpb;
850 int16_t tag;
851 uint8_t i,ie,iface,iea,ieb;
852
853 iel = lists[l] / 4;
854 iface = lists[l] % 4;
855 pt = &mesh->tetra[iel];
856 iea = ieb = 0;
857
858 assert ( pt->xt && "tetra with boundary face has a xtetra");
859
860 /* For each bdy face that contains ip0, store the index of the 2 edges
861 * passing through \a ip0 in \a iea and \a ieb. */
862 for (i=0; i<3; i++) {
863 ie = MMG5_iarf[iface][i]; //index in tet of edge i on face iface
864 if ( (pt->v[MMG5_iare[ie][0]] == ip0) || (pt->v[MMG5_iare[ie][1]] == ip0) ) {
865 if ( !iea )
866 iea = ie;
867 else
868 ieb = ie;
869 }
870 }
871 /* In current face (\a iface), store in \a iptmpa the global index of the
872 * second vertex of edge \a iea (first vertex being \a ip0). */
873 if ( pt->v[MMG5_iare[iea][0]] != ip0 )
874 iptmpa = pt->v[MMG5_iare[iea][0]];
875 else {
876 assert(pt->v[MMG5_iare[iea][1]] != ip0);
877 iptmpa = pt->v[MMG5_iare[iea][1]];
878 }
879 /* In current face (\a iface), store in \a iptmpb the global index of the
880 * second vertex of edge \a ieb (first vertex being \a ip0). */
881 if ( pt->v[MMG5_iare[ieb][0]] != ip0 )
882 iptmpb = pt->v[MMG5_iare[ieb][0]];
883 else {
884 assert(pt->v[MMG5_iare[ieb][1]] != ip0);
885 iptmpb = pt->v[MMG5_iare[ieb][1]];
886 }
887
888 /* Search if the edge ip0-iptmpa is the edge at the interface with previous
889 * triangle. */
890 if ( (iptmpa == *ipa) || (iptmpa == *ipb) ) {
891 tag = mesh->xtetra[pt->xt].tag[iea];
892 if ( edgTag & tag ) {
893 /* The featured edge has been found. End of ball processing. */
894 *ip = iptmpa;
895 return 1;
896 }
897 }
898
899 /* Search if the edge ip0-iptmpb is the edge at the interface with previous
900 * triangle. */
901 if ( (iptmpb == *ipa) || (iptmpb == *ipb) ) {
902 tag = mesh->xtetra[pt->xt].tag[ieb];
903 if ( edgTag & tag ) {
904 /* The featured edge has been found. End of ball processing. */
905 *ip = iptmpb;
906 return 1;
907 }
908 }
909
910 /* Update face vertices for next item processing */
911 *ipa = iptmpa;
912 *ipb = iptmpb;
913
914 return 0;
915}
916
936int MMG3D_curveEndingPts(MMG5_pMesh mesh,MMG5_int *lists,int ilists,
937 const int16_t edgTag, MMG5_int ip0,MMG5_int *ip1,
938 MMG5_int *ip2) {
939 MMG5_pTetra pt;
940 MMG5_int iel,ipa,ipb;
941 int l;
942 uint8_t i,iface;
943
949 /* Get first edge to initialize the loop: \a ip0 is the global idx of the
950 * point we want to move, store in \a ipa and \a ipb the global indices of the
951 * 2 other vertices of the boundary face from which we start. When processing
952 * next triangle we will find either the edge ip0-ipa, or ip0-ipb, this will
953 * be the first edge that we will check. */
954 iel = lists[0]/4;
955 iface = lists[0]%4;
956 pt = &mesh->tetra[iel];
957 ipa = ipb = 0;
958 for (i=0; i<3; i++) {
959 if ( pt->v[MMG5_idir[iface][i]] != ip0 ) {
960 if ( !ipa )
961 ipa = pt->v[MMG5_idir[iface][i]];
962 else
963 ipb = pt->v[MMG5_idir[iface][i]];
964 }
965 }
966 assert(ipa && ipb);
967
968 /* Travel surfacic list of \a ip0 and search if the edge at interface of
969 * boundary triangles stored in lists[l] and lists[l-1] belongs to our curve
970 * (\a edgTag edge). */
971 for (l=1; l<ilists; l++) {
972 if ( MMG3D_curveEndingPts_chkEdg(mesh,lists,l,ip0,&ipa,&ipb,edgTag,ip1) ) {
973 break;
974 }
975 }
976
983 /* Get first edge to initialize the loop: \a ip0 is the global idx of the
984 * point we want to move, store in \a ipa and \a ipb the global indices of the
985 * 2 other vertices of the boundary face from which we start. When processing
986 * next triangle we will find either the edge ip0-ipa, or ip0-ipb, this will
987 * be the first edge that we will check. */
988 iel = lists[0]/4;
989 iface = lists[0]%4;
990 pt = &mesh->tetra[iel];
991 ipa = ipb = 0;
992 for (i=0; i<3; i++) {
993 if ( pt->v[MMG5_idir[iface][i]] != ip0 ) {
994 if ( !ipa )
995 ipa = pt->v[MMG5_idir[iface][i]];
996 else
997 ipb = pt->v[MMG5_idir[iface][i]];
998 }
999 }
1000 assert(ipa && ipb);
1001
1002 /* Travel surfacic list of \a ip0 and search if the edge at interface of
1003 * boundary triangles stored in lists[l] and lists[l+1] belongs to our curve
1004 * (\a edgTag edge). */
1005 for (l=ilists-1; l>0; l--) {
1006 if ( MMG3D_curveEndingPts_chkEdg(mesh,lists,l,ip0,&ipa,&ipb,edgTag,ip2) ) {
1007 break;
1008 }
1009 }
1010
1011 /* Check that we have found two distinct ending points */
1012 if ( !(*ip1) ) {
1013 return 0;
1014 }
1015 if ( !(*ip2) ) {
1016 return 0;
1017 }
1018 if ( (*ip1) == (*ip2) ) {
1019 return 0;
1020 }
1021
1022 return 1;
1023}
1024
1049 MMG3D_pPROctree PROctree, int64_t *listv,
1050 int ilistv,int improve,MMG5_pPoint p0,
1051 MMG5_int ip0,uint8_t isrid,double o[3],
1052 double no[3],double no2[3],double to[3]) {
1053
1054 MMG5_pTetra pt,pt0;
1055 MMG5_pxPoint pxp;
1056 double calold,calnew,callist[MMG3D_LMAX+2];
1057 MMG5_int iel;
1058 int l;
1059 int8_t i0;
1060
1062 calold = calnew = DBL_MAX;
1063
1064 for( l=0 ; l<ilistv ; l++ ){
1065 iel = listv[l] / 4;
1066 i0 = listv[l] % 4;
1067 pt = &mesh->tetra[iel];
1068 pt0 = &mesh->tetra[0];
1069 memcpy(pt0,pt,sizeof(MMG5_Tetra));
1070 pt0->v[i0] = 0;
1071 calold = MG_MIN(calold, pt->qual);
1072 callist[l] = MMG5_orcal(mesh,met,0);
1073 if (callist[l] < MMG5_NULKAL) {
1074 return 0;
1075 }
1076 calnew = MG_MIN(calnew,callist[l]);
1077 }
1078 if ((calold < MMG5_EPSOK && calnew <= calold) ||
1079 (calnew < MMG5_EPSOK) || (calnew <= 0.3*calold)) {
1080 return 0;
1081 } else if (improve && calnew < calold) {
1082 return 0;
1083 }
1084
1086 if ( PROctree ) {
1087 MMG3D_movePROctree(mesh, PROctree, ip0, o, p0->c);
1088 }
1089
1090 p0->c[0] = o[0];
1091 p0->c[1] = o[1];
1092 p0->c[2] = o[2];
1093
1094 pxp = &mesh->xpoint[p0->xp];
1095 pxp->n1[0] = no[0];
1096 pxp->n1[1] = no[1];
1097 pxp->n1[2] = no[2];
1098
1099 p0->n[0] = to[0];
1100 p0->n[1] = to[1];
1101 p0->n[2] = to[2];
1102
1103 if ( isrid ) {
1104 /* Copy the second normal for ridge point */
1105 pxp->n2[0] = no2[0];
1106 pxp->n2[1] = no2[1];
1107 pxp->n2[2] = no2[2];
1108 }
1109
1110 for( l=0 ; l<ilistv ; l++ ){
1111 (&mesh->tetra[listv[l]/4])->qual = callist[l];
1112 (&mesh->tetra[listv[l]/4])->mark = mesh->mark;
1113 }
1114 return 1;
1115}
1116
1143 MMG5_int ip1,MMG5_int ip2,double ll1old,double ll2old,
1144 uint8_t isrid,const double step,
1145 double o[3],double no[3],
1146 double no2[3],double to[3],
1147 const int16_t edgTag) {
1148
1149 MMG5_int ip;
1150
1152 if ( ll1old < ll2old ) {
1153 /* move towards p2 */
1154 ip = ip2;
1155 }
1156 else if ( ll1old > ll2old ) {
1157 /* move towards p1 */
1158 ip = ip1;
1159 }
1160 else {
1161 return 0;
1162 }
1163
1165 if ( MG_NOM & edgTag ) {
1166 if ( !(MMG5_BezierNom(mesh,ip0,ip,step,o,no,to)) ) {
1167 return 0;
1168 }
1169 }
1170 else if ( MG_GEO & edgTag ) {
1171 // Remark: Singular points are required so following assertion should be
1172 // verified in the entire function. Keep the test here to make easier
1173 // debugging/understanding when passing here.
1174 assert ( (!MG_SIN(mesh->point[ip0].tag)) &&
1175 "BezierRidge don't work if both ip0 and ip are singular" );
1176 if ( !(MMG5_BezierRidge(mesh,ip0,ip,step,o,no,no2,to)) ) {
1177 return 0;
1178 }
1179 }
1180 else if ( MG_REF & edgTag ) {
1181 if ( !(MMG5_BezierRef(mesh,ip0,ip,step,o,no,to)) ) {
1182 return 0;
1183 }
1184 }
1185 else {
1186 assert ( 0 && "Unexpected edge tag in this function");
1187 return 0;
1188 }
1189
1191 MMG5_pPoint ppt0 = &mesh->point[0];
1192 ppt0->c[0] = o[0];
1193 ppt0->c[1] = o[1];
1194 ppt0->c[2] = o[2];
1195 ppt0->tag = p0->tag;
1196 ppt0->ref = p0->ref;
1197
1198
1199 MMG5_int nxp = mesh->xp + 1;
1200 if ( nxp > mesh->xpmax ) {
1202 "larger xpoint table",
1203 return 0);
1204 }
1205 ppt0->xp = nxp;
1206 MMG5_pxPoint pxp = &mesh->xpoint[nxp];
1207 memcpy(pxp,&(mesh->xpoint[p0->xp]),sizeof(MMG5_xPoint));
1208
1209 ppt0->n[0] = to[0];
1210 ppt0->n[1] = to[1];
1211 ppt0->n[2] = to[2];
1212
1213 pxp->n1[0] = no[0];
1214 pxp->n1[1] = no[1];
1215 pxp->n1[2] = no[2];
1216
1217 if ( isrid ) {
1218 /* Copy the second normal for ridge point */
1219 pxp->n2[0] = no2[0];
1220 pxp->n2[1] = no2[1];
1221 pxp->n2[2] = no2[2];
1222 }
1223
1224 return ip;
1225}
1226
1227
1254static inline
1256 int ilistv, MMG5_int *lists, int ilists,int improve,const int16_t edgTag){
1257 MMG5_pTetra pt;
1258 MMG5_pxTetra pxt;
1259 MMG5_pPoint p0,p1,p2;
1260 MMG5_Tria tt;
1261 MMG5_pPar par;
1262 double ll1old,ll2old,o[3],no[3],no2[3],to[3];
1263 double calold,calnew,caltmp,hmax,hausd;
1264 MMG5_int iel,ip0,ip1,ip2,ip;
1265 int l;
1266 int isloc,j;
1267 uint8_t i,iface,isrid;
1268
1269 ip1 = ip2 = 0;
1270 pt = &mesh->tetra[listv[0]/4];
1271 ip0 = pt->v[listv[0]%4];
1272 p0 = &mesh->point[ip0];
1273
1276 isrid = ((MG_GEO & edgTag) && !(MG_NOM & edgTag));
1277
1278 assert ( edgTag & p0->tag );
1279
1282 int ier = MMG3D_curveEndingPts(mesh,lists,ilists,edgTag,ip0,&ip1,&ip2);
1283 if ( !ier ) {
1284 return 0;
1285 }
1286
1293 p1 = &mesh->point[ip1];
1294 p2 = &mesh->point[ip2];
1295
1296 ll1old = (p1->c[0] -p0->c[0])* (p1->c[0] -p0->c[0]) \
1297 + (p1->c[1] -p0->c[1])* (p1->c[1] -p0->c[1]) \
1298 + (p1->c[2] -p0->c[2])* (p1->c[2] -p0->c[2]);
1299 ll2old = (p2->c[0] -p0->c[0])* (p2->c[0] -p0->c[0]) \
1300 + (p2->c[1] -p0->c[1])* (p2->c[1] -p0->c[1]) \
1301 + (p2->c[2] -p0->c[2])* (p2->c[2] -p0->c[2]);
1302
1305 ip = MMG3D_movbdycurvept_newPosForSimu( mesh,p0,ip0,ip1,ip2,ll1old,ll2old,
1306 isrid,MMG3D_MOVSTEP,o,no,no2,to,edgTag );
1307 if ( !ip ) {
1308 return 0;
1309 }
1310
1320 calold = calnew = DBL_MAX;
1321 for( l=0 ; l<ilists ; l++ ){
1322 iel = lists[l] / 4;
1323 iface = lists[l] % 4;
1324
1325 assert( 0<=iface && iface<4 && "unexpected local face idx");
1326 MMG5_tet2tri(mesh,iel,iface,&tt);
1327 caltmp = MMG5_caltri(mesh,met,&tt);
1328 calold = MG_MIN(calold,caltmp);
1329
1330 for( i=0 ; i<3 ; i++ ) {
1331 if ( tt.v[i] == ip0 ) {
1332 break;
1333 }
1334 }
1335 assert(i<3);
1336 if ( i==3 ) {
1337 return 0;
1338 }
1339
1340 tt.v[i] = 0;
1341
1342 caltmp = MMG5_caltri(mesh,met,&tt);
1343 if ( caltmp < MMG5_EPSD2 ) {
1344 /* We don't check the input triangle qualities, thus we may have a very
1345 * bad triangle in our mesh */
1346 return 0;
1347 }
1348 calnew = MG_MIN(calnew,caltmp);
1349
1350 /* Local parameters for tt and iel */
1351 pt = &mesh->tetra[iel];
1352 pxt = &mesh->xtetra[pt->xt];
1353
1354 hmax = mesh->info.hmax;
1355 hausd = mesh->info.hausd;
1356
1357 isloc = 0;
1358 if ( mesh->info.parTyp & MG_Tetra ) {
1359 for ( j=0; j<mesh->info.npar; ++j ) {
1360 par = &mesh->info.par[j];
1361
1362 if ( par->elt != MMG5_Tetrahedron ) continue;
1363 if ( par->ref != pt->ref ) continue;
1364
1365 hmax = par->hmax;
1366 hausd = par->hausd;
1367 isloc = 1;
1368 break;
1369 }
1370 }
1371 if ( mesh->info.parTyp & MG_Tria ) {
1372 if ( isloc ) {
1373 for ( j=0; j<mesh->info.npar; ++j ) {
1374 par = &mesh->info.par[j];
1375
1376 if ( par->elt != MMG5_Triangle ) continue;
1377 if ( par->ref != tt.ref ) continue;
1378
1379 hmax = MG_MIN(hmax,par->hmax);
1380 hausd = MG_MIN(hausd,par->hausd);
1381 break;
1382 }
1383 }
1384 else {
1385 for ( j=0; j<mesh->info.npar; ++j ) {
1386 par = &mesh->info.par[j];
1387
1388 if ( par->elt != MMG5_Triangle ) continue;
1389 if ( par->ref != tt.ref ) continue;
1390
1391 hmax = par->hmax;
1392 hausd = par->hausd;
1393 isloc = 1;
1394 break;
1395 }
1396 }
1397 }
1398
1399 if ( MMG5_chkedg(mesh,&tt,MG_GET(pxt->ori,iface),hmax,hausd,isloc) > 0 ) {
1400 memset(&mesh->xpoint[mesh->point[0].xp],0,sizeof(MMG5_xPoint));
1401 return 0;
1402 }
1403 }
1404 if ( calold < MMG5_EPSOK && calnew <= calold ) {
1405 return 0;
1406 }
1407 else if ( calnew < calold ) {
1408 return 0;
1409 }
1410 memset(&mesh->xpoint[mesh->point[0].xp],0,sizeof(MMG5_xPoint));
1411
1415 ier = MMG3D_movbdycurvept_chckAndUpdate(mesh,met,PROctree,listv,ilistv,
1416 improve,p0,ip0,isrid,o,no,no2,to);
1417
1418 return ier;
1419}
1420
1439 int ilistv, MMG5_int *lists, int ilists,int improve){
1440
1441 return MMG3D_movbdycurvept_iso(mesh,met,PROctree,listv,ilistv,lists,ilists,improve,MG_REF);
1442}
1443
1463 int ilistv, MMG5_int *lists, int ilists,int improve){
1464
1465 return MMG3D_movbdycurvept_iso(mesh,met,PROctree,listv,ilistv,lists,ilists,improve,MG_NOM);
1466}
1467
1483 int ilistv, int improve){
1484 MMG5_pTetra pt,pt0;
1485 MMG5_pxTetra pxt;
1486 MMG5_pPoint p0,p1,p2,ppt0;
1487 double step,ll1old,ll2old,calold,calnew,callist[MMG3D_LMAX+2];
1488 double o[3],no[3],to[3];
1489 int l;
1490 MMG5_int ip0,ip1,ip2,ip,iel,ipa;
1491 int8_t i,i0,ie;
1492
1493 step = 0.1;
1494 ip1 = ip2 = 0;
1495 pt = &mesh->tetra[listv[0]/4];
1496 ip0 = pt->v[listv[0]%4];
1497 p0 = &mesh->point[ip0];
1498
1499 assert ( (p0->tag & MG_NOM) && p0->xp && mesh->xpoint[p0->xp].nnor );
1500
1501 /* Recover the two ending points of the underlying non manifold curve */
1502 for (l=0; l<ilistv; l++) {
1503 iel = listv[l] / 4;
1504 i0 = listv[l] % 4;
1505 pt = &mesh->tetra[iel];
1506 if ( !pt->xt ) continue;
1507 pxt = &mesh->xtetra[pt->xt];
1508
1509 for (i=0; i<3; i++) {
1510 ie = MMG5_arpt[i0][i];
1511 if ( pxt->tag[ie] & MG_NOM ) {
1512 ipa = ( ip0 == pt->v[MMG5_iare[ie][0]] ) ? pt->v[MMG5_iare[ie][1]] : pt->v[MMG5_iare[ie][0]];
1513 if ( !ip1 ) ip1 = ipa;
1514 else if ( !ip2 && ipa != ip1 ) ip2 = ipa;
1515 }
1516 }
1517 }
1518 if ( !(ip1 && ip2 && (ip1 != ip2)) ) return 0;
1519
1520 /* At this point, we get the point extremities ip1, ip2 of the non manifold curve passing through ip0 */
1521 p1 = &mesh->point[ip1];
1522 p2 = &mesh->point[ip2];
1523
1524 ll1old = (p1->c[0] -p0->c[0])* (p1->c[0] -p0->c[0]) \
1525 + (p1->c[1] -p0->c[1])* (p1->c[1] -p0->c[1]) \
1526 + (p1->c[2] -p0->c[2])* (p1->c[2] -p0->c[2]);
1527 ll2old = (p2->c[0] -p0->c[0])* (p2->c[0] -p0->c[0]) \
1528 + (p2->c[1] -p0->c[1])* (p2->c[1] -p0->c[1]) \
1529 + (p2->c[2] -p0->c[2])* (p2->c[2] -p0->c[2]);
1530
1531 if ( ll1old < ll2old ) { //move towards p2
1532 ip = ip2;
1533 }
1534 else {
1535 ip = ip1;
1536 }
1537
1538 /* Compute support of the associated edge, and features of the new position */
1539 if ( !(MMG5_BezierNom(mesh,ip0,ip,step,o,no,to)) ) return 0;
1540
1541 /* Test : check whether all volumes remain positive with new position of the point */
1542 ppt0 = &mesh->point[0];
1543 ppt0->c[0] = o[0];
1544 ppt0->c[1] = o[1];
1545 ppt0->c[2] = o[2];
1546 ppt0->tag = p0->tag;
1547 ppt0->ref = p0->ref;
1548
1549 calold = calnew = DBL_MAX;
1550 for( l=0 ; l<ilistv ; l++ ){
1551 iel = listv[l] / 4;
1552 i0 = listv[l] % 4;
1553 pt = &mesh->tetra[iel];
1554 pt0 = &mesh->tetra[0];
1555 memcpy(pt0,pt,sizeof(MMG5_Tetra));
1556 pt0->v[i0] = 0;
1557 calold = MG_MIN(calold, pt->qual);
1558 callist[l]= MMG5_orcal(mesh,met,0);
1559 if (callist[l] < MMG5_NULKAL) {
1560 return 0;
1561 }
1562 calnew = MG_MIN(calnew,callist[l]);
1563 }
1564 if ((calold < MMG5_EPSOK && calnew <= calold) ||
1565 (calnew < MMG5_EPSOK) || (calnew <= 0.3*calold)) {
1566 return 0;
1567 } else if (improve && calnew < calold) {
1568 return 0;
1569 }
1570
1571 /* Update coordinates and tangent for new point */
1572 if ( PROctree )
1573 MMG3D_movePROctree(mesh, PROctree, ip0, o, p0->c);
1574
1575 p0->c[0] = o[0];
1576 p0->c[1] = o[1];
1577 p0->c[2] = o[2];
1578
1579 p0->n[0] = to[0];
1580 p0->n[1] = to[1];
1581 p0->n[2] = to[2];
1582
1583 for(l=0; l<ilistv; l++){
1584 (&mesh->tetra[listv[l]/4])->qual = callist[l];
1585 (&mesh->tetra[listv[l]/4])->mark = mesh->mark;
1586 }
1587
1588 return 1;
1589}
1590
1607 int ilistv,MMG5_int *lists,int ilists,int improve) {
1608
1609 return MMG3D_movbdycurvept_iso(mesh,met,PROctree,listv,ilistv,lists,ilists,improve,MG_GEO);
1610}
1611
1627 MMG5_pTetra pt,pt1;
1628 MMG5_pPoint ppa,ppb,p1,p2,p3;
1629 int j,iter,maxiter,l,lon;
1630 int64_t list[MMG3D_LMAX+2];
1631 MMG5_int ipb,iadr,i1,i2,i3,iel;
1632 double *mp,coe,qualtet[MMG3D_LMAX+2];
1633 double ax,ay,az,bx,by,bz,nx,ny,nz,dd,len,qual,oldc[3];
1634
1635 assert(k);
1636 assert(ib<4);
1637 pt = &mesh->tetra[k];
1638 ppa = &mesh->point[pt->v[ib]];
1639
1640 if ( (ppa->tag & MG_BDY) || (ppa->tag & MG_REQ) ) {
1641 return 0;
1642 }
1643
1644 iadr = pt->v[ib]*sol->size + 0;
1645 mp = &sol->m[iadr];
1646
1647 /*compute normal*/
1648 i1 = pt->v[MMG5_idir[ib][0]];
1649 i2 = pt->v[MMG5_idir[ib][1]];
1650 i3 = pt->v[MMG5_idir[ib][2]];
1651 p1 = &mesh->point[i1];
1652 p2 = &mesh->point[i2];
1653 p3 = &mesh->point[i3];
1654
1655 ax = p3->c[0] - p1->c[0];
1656 ay = p3->c[1] - p1->c[1];
1657 az = p3->c[2] - p1->c[2];
1658
1659 bx = p2->c[0] - p1->c[0];
1660 by = p2->c[1] - p1->c[1];
1661 bz = p2->c[2] - p1->c[2];
1662
1663 nx = (ay*bz - az*by);
1664 ny = (az*bx - ax*bz);
1665 nz = (ax*by - ay*bx);
1666
1667 dd = sqrt(nx*nx+ny*ny+nz*nz);
1668 dd = 1./dd;
1669 nx *=dd;
1670 ny *=dd;
1671 nz *=dd;
1672 len = 0;
1673 for (j=0; j<3; j++) {
1674 ipb = pt->v[ MMG5_idir[ib][j] ];
1675 ppb = &mesh->point[ipb];
1676
1677 ax = ppb->c[0] - ppa->c[0];
1678 ay = ppb->c[1] - ppa->c[1];
1679 az = ppb->c[2] - ppa->c[2];
1680
1681 dd = mp[0]*ax*ax + mp[3]*ay*ay + mp[5]*az*az \
1682 + 2.0*(mp[1]*ax*ay + mp[2]*ax*az + mp[4]*ay*az);
1683 assert(dd>0);
1684 len += sqrt(dd);
1685 }
1686
1687 dd = 1.0 / (double)3.;
1688 len = 1.0 / len;
1689 len *= dd;
1690 memcpy(oldc,ppa->c,3*sizeof(double));
1691
1692 lon = MMG5_boulevolp(mesh,k,ib,&list[0]);
1693
1694 coe = 1.;
1695 iter = 0;
1696 maxiter = 20;
1697 do {
1698 ppa->c[0] = oldc[0] + coe * nx * len;
1699 ppa->c[1] = oldc[1] + coe * ny * len;
1700 ppa->c[2] = oldc[2] + coe * nz * len;
1701
1702 for (l=0; l<lon; l++) {
1703 iel = list[l] / 4 ;
1704 pt1 = &mesh->tetra[iel];
1705
1706 qual = MMG5_caltet(mesh,sol,pt1);
1707 /*warning if we increase the coefficient (ex 1.4), the mesh quality becomes poor very quickly*/
1708 if ( qual*1.01 <= pt1->qual) break;
1709 qualtet[l] = qual;
1710
1711 }
1712 if ( l >= lon ) break;
1713 coe *= 0.5;
1714 }
1715 while ( ++iter <= maxiter );
1716 if ( iter > maxiter) {
1717 memcpy(ppa->c,oldc,3*sizeof(double));
1718 return 0;
1719 }
1720
1721 for (l=0; l<lon; l++) {
1722 iel = list[l] / 4;
1723 pt1 = &mesh->tetra[iel];
1724 pt1->qual = qualtet[l];
1725 pt1->mark = mesh->mark;
1726 }
1727 return 1;
1728
1729}
1730
1747 MMG5_pTetra pt,pt1;
1748 MMG5_pPoint ppa,ppb,p1,p2,p3;
1749 int j,iter,maxiter,l,lon;
1750 int64_t list[MMG3D_LMAX+2];
1751 MMG5_int ipb,iel,i1,i2,i3;
1752 double coe,crit,qualtet[MMG3D_LMAX+2];
1753 double ax,ay,az,bx,by,bz,nx,ny,nz,dd,len,qual,oldc[3],oldp[3];
1754
1755 assert(k);
1756 assert(ib<4);
1757 pt = &mesh->tetra[k];
1758
1759 ppa = &mesh->point[pt->v[ib]];
1760
1761 if ( ppa->tag & MG_BDY || (ppa->tag & MG_REQ) ) {
1762 return 0;
1763 }
1764
1765 /*compute normal*/
1766 i1 = pt->v[MMG5_idir[ib][0]];
1767 i2 = pt->v[MMG5_idir[ib][1]];
1768 i3 = pt->v[MMG5_idir[ib][2]];
1769 p1 = &mesh->point[i1];
1770 p2 = &mesh->point[i2];
1771 p3 = &mesh->point[i3];
1772
1773 ax = p3->c[0] - p1->c[0];
1774 ay = p3->c[1] - p1->c[1];
1775 az = p3->c[2] - p1->c[2];
1776
1777 bx = p2->c[0] - p1->c[0];
1778 by = p2->c[1] - p1->c[1];
1779 bz = p2->c[2] - p1->c[2];
1780
1781 nx = (ay*bz - az*by);
1782 ny = (az*bx - ax*bz);
1783 nz = (ax*by - ay*bx);
1784
1785 dd = sqrt(nx*nx+ny*ny+nz*nz);
1786 dd = 1./dd;
1787 nx *=dd;
1788 ny *=dd;
1789 nz *=dd;
1790 len = 0;
1791 for (j=0; j<3; j++) {
1792 ipb = pt->v[ MMG5_idir[ib][j] ];
1793 ppb = &mesh->point[ipb];
1794
1795 ax = ppb->c[0] - ppa->c[0];
1796 ay = ppb->c[1] - ppa->c[1];
1797 az = ppb->c[2] - ppa->c[2];
1798
1799 dd = sqrt(ax*ax +ay*ay +az*az);
1800 len += dd;
1801 }
1802
1803 dd = 1.0 / (double)3.;
1804 len *= dd;
1805
1806
1807 memcpy(oldp,ppa->c,3*sizeof(double));
1808 oldc[0] = 1./3.*(p1->c[0]+p2->c[0]+p3->c[0]);
1809 oldc[1] = 1./3.*(p1->c[1]+p2->c[1]+p3->c[1]);
1810 oldc[2] = 1./3.*(p1->c[2]+p2->c[2]+p3->c[2]);
1811
1812 lon = MMG5_boulevolp(mesh,k,ib,&list[0]);
1813
1814 if ( !lon ) return 0;
1815
1816 /*vol crit*/
1817 crit = MMG5_orvol(mesh->point,pt->v);
1818 for (l=1; l<lon; l++) {
1819 iel = list[l] / 4;
1820 pt1 = &mesh->tetra[iel];
1821 if ( pt1->qual < crit) crit = pt1->qual;
1822 }
1823 coe = 0.471404;//2.12132; //3/sqrt(2) : hauteur d'un tetra reg de cote c : c*sqrt(2)/3
1824 iter = 0;
1825 maxiter = 10;
1826 do {
1827 ppa->c[0] = oldc[0] + coe * nx * len;
1828 ppa->c[1] = oldc[1] + coe * ny * len;
1829 ppa->c[2] = oldc[2] + coe * nz * len;
1830 for (l=0; l<lon; l++) {
1831 iel = list[l] / 4;
1832 pt1 = &mesh->tetra[iel];
1833 qual = MMG5_caltet(mesh,sol,pt1);
1834 if ( qual < crit ) {
1835 break;
1836 }
1837 qualtet[l] = qual;
1838
1839 }
1840 if ( l >= lon ) break;
1841 coe *= 0.5;
1842 }
1843 while ( ++iter <= maxiter );
1844 if ( iter > maxiter) {
1845 memcpy(ppa->c,oldp,3*sizeof(double));
1846 return 0;
1847 }
1848
1849 for (l=0; l<lon; l++) {
1850 iel = list[l] / 4;
1851 pt1 = &mesh->tetra[iel];
1852 pt1->qual = qualtet[l];
1853 pt1->mark = mesh->mark;
1854 }
1855 return 1;
1856
1857}
1858
1876 MMG5_pTetra pt,pt1;
1877 MMG5_pPoint ppa,ppb,p1,p2,p3;
1878 int j,iter,maxiter,l,lon;
1879 int64_t list[MMG3D_LMAX+2];
1880 MMG5_int ipb,iadr,iel,i1,i2,i3;
1881 double hp,coe,crit,qualtet[MMG3D_LMAX+2];;
1882 double ax,ay,az,bx,by,bz,nx,ny,nz,dd,len,qual,oldc[3];
1883
1884 assert(k);
1885 assert(ib<4);
1886 pt = &mesh->tetra[k];
1887
1888 ppa = &mesh->point[pt->v[ib]];
1889 if ( (ppa->tag & MG_BDY) || (ppa->tag & MG_REQ) ) {
1890 return 0;
1891 }
1892
1893 iadr = (pt->v[ib])*sol->size;
1894 hp = sol->m[iadr];
1895
1896 /*compute normal*/
1897 i1 = pt->v[MMG5_idir[ib][0]];
1898 i2 = pt->v[MMG5_idir[ib][1]];
1899 i3 = pt->v[MMG5_idir[ib][2]];
1900 p1 = &mesh->point[i1];
1901 p2 = &mesh->point[i2];
1902 p3 = &mesh->point[i3];
1903
1904 ax = p3->c[0] - p1->c[0];
1905 ay = p3->c[1] - p1->c[1];
1906 az = p3->c[2] - p1->c[2];
1907
1908 bx = p2->c[0] - p1->c[0];
1909 by = p2->c[1] - p1->c[1];
1910 bz = p2->c[2] - p1->c[2];
1911
1912 nx = (ay*bz - az*by);
1913 ny = (az*bx - ax*bz);
1914 nz = (ax*by - ay*bx);
1915
1916 dd = sqrt(nx*nx+ny*ny+nz*nz);
1917 dd = 1./dd;
1918 nx *=dd;
1919 ny *=dd;
1920 nz *=dd;
1921 len = 0;
1922 for (j=0; j<3; j++) {
1923 ipb = pt->v[ MMG5_idir[ib][j] ];
1924 ppb = &mesh->point[ipb];
1925
1926 ax = ppb->c[0] - ppa->c[0];
1927 ay = ppb->c[1] - ppa->c[1];
1928 az = ppb->c[2] - ppa->c[2];
1929
1930 dd = sqrt(ax*ax +ay*ay +az*az);
1931 len += dd/hp;
1932 }
1933
1934 dd = 1.0 / (double)3.;
1935 len *= dd;
1936 if(len > 0.) len = 1.0 / len;
1937
1938 memcpy(oldc,ppa->c,3*sizeof(double));
1939
1940 lon = MMG5_boulevolp(mesh,k,ib,&list[0]);
1941
1942 if(!lon) return 0;
1943
1944 /*qual crit*/
1945 crit = pt->qual;
1946 for (l=1; l<lon; l++) {
1947 iel = list[l] / 4;
1948 pt1 = &mesh->tetra[iel];
1949 if ( pt1->qual < crit )
1950 crit = pt1->qual;
1951 }
1952
1953 crit *= 1.01;
1954 coe = 1.;
1955 iter = 0;
1956 maxiter = 20;
1957 do {
1958
1959 ppa->c[0] = oldc[0] + coe * nx * len;
1960 ppa->c[1] = oldc[1] + coe * ny * len;
1961 ppa->c[2] = oldc[2] + coe * nz * len;
1962 for (l=0; l<lon; l++) {
1963 iel = list[l] / 4;
1964 pt1 = &mesh->tetra[iel];
1965 qual = MMG5_caltet(mesh,sol,pt1);
1966 if ( qual < crit ) break;
1967 qualtet[l] = qual;
1968
1969 }
1970 if ( l >= lon ) break;
1971 coe *= 0.5;
1972 }
1973 while ( ++iter <= maxiter );
1974 if ( iter > maxiter) {
1975 memcpy(ppa->c,oldc,3*sizeof(double));
1976 return 0;
1977 }
1978
1979 for (l=0; l<lon; l++) {
1980 iel = list[l] / 4;
1981 pt1 = &mesh->tetra[iel];
1982 pt1->qual = qualtet[l];
1983 pt1->mark = mesh->mark;
1984 }
1985 return 1;
1986
1987}
int ier
MMG5_pMesh MMG5_pSol * sol
if(!ier) exit(EXIT_FAILURE)
MMG5_pMesh * mesh
int MMG3D_movePROctree(MMG5_pMesh mesh, MMG3D_pPROctree q, MMG5_int no, double *newVer, double *oldVer)
Definition: PRoctree_3d.c:225
int MMG3D_bezierInt(MMG5_pBezier pb, double uv[2], double o[3], double no[3], double to[3])
Definition: bezier_3d.c:608
int MMG5_boulevolp(MMG5_pMesh mesh, MMG5_int start, int ip, int64_t *list)
Definition: boulep_3d.c:54
static double MMG5_lenedg_iso(MMG5_pMesh mesh, MMG5_pSol met, int ia, MMG5_pTetra pt)
static double MMG5_orcal(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int iel)
API headers for the mmg3d library.
#define MMG3D_LMAX
Definition: libmmg3d.h:58
#define MMG3D_MOVSTEP
#define MMG3D_LOPTS
static const uint8_t MMG5_arpt[4][3]
arpt[i]: edges passing through vertex i
int MMG5_BezierRidge(MMG5_pMesh mesh, MMG5_int ip0, MMG5_int ip1, double s, double *o, double *no1, double *no2, double *to)
Definition: tools_3d.c:145
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
static const uint8_t MMG5_iare[6][2]
vertices of extremities of the edges of the tetra
int MMG5_BezierRef(MMG5_pMesh mesh, MMG5_int ip0, MMG5_int ip1, double s, double *o, double *no, double *to)
Definition: tools_3d.c:352
void MMG5_tet2tri(MMG5_pMesh mesh, MMG5_int k, int8_t ie, MMG5_Tria *ptt)
Definition: mmg3d1.c:102
#define MMG3D_LOPTL
int MMG5_BezierNom(MMG5_pMesh mesh, MMG5_int ip0, MMG5_int ip1, double s, double *o, double *no, double *to)
Definition: tools_3d.c:527
static const uint8_t MMG5_idir[4][3]
idir[i]: vertices of face opposite to vertex i
@ MMG5_Tetrahedron
Definition: libmmgtypes.h:227
@ MMG5_Triangle
Definition: libmmgtypes.h:226
#define MG_Tria
#define MG_REQ
#define MG_GEO
#define MMG5_EPSOK
#define MMG5_SQR32
#define MG_EOK(pt)
#define MG_MIN(a, b)
#define MMG5_GAP
#define MG_EDG(tag)
double MMG5_det3pt1vec(double c0[3], double c1[3], double c2[3], double v[3])
Definition: tools.c:920
double MMG5_det4pt(double c0[3], double c1[3], double c2[3], double c3[3])
Definition: tools.c:932
#define MMG5_TAB_RECALLOC(mesh, ptr, initSize, wantedGap, type, message, law)
#define MG_SIN(tag)
#define MMG5_NULKAL
#define MG_Tetra
#define MG_GET(flag, bit)
double MMG5_orvol(MMG5_pPoint point, MMG5_int *v)
Definition: tools.c:951
int MMG5_rotmatrix(double n[3], double r[3][3])
Definition: tools.c:467
#define MG_BDY
#define MG_NOM
#define MMG5_EPSD2
#define MG_REF
int MMG5_movbdynomintpt_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *listv, int ilistv, int improve)
Definition: movpt_3d.c:1482
int MMG5_movbdyridpt_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *listv, int ilistv, MMG5_int *lists, int ilists, int improve)
Definition: movpt_3d.c:1606
int MMG3D_movnormal_iso(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int k, int ib)
Definition: movpt_3d.c:1746
int MMG3D_movbdycurvept_chckAndUpdate(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *listv, int ilistv, int improve, MMG5_pPoint p0, MMG5_int ip0, uint8_t isrid, double o[3], double no[3], double no2[3], double to[3])
Definition: movpt_3d.c:1048
int MMG3D_rotate_surfacicBall(MMG5_pMesh mesh, MMG5_int *lists, int ilists, MMG5_int ip0, double r[3][3], double *lispoi)
Definition: movpt_3d.c:350
int MMG3D_movv_ani(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int k, int ib)
Definition: movpt_3d.c:1626
static int MMG3D_curveEndingPts_chkEdg(MMG5_pMesh mesh, MMG5_int *lists, int l, MMG5_int ip0, MMG5_int *ipa, MMG5_int *ipb, const int16_t edgTag, MMG5_int *ip)
Definition: movpt_3d.c:845
int MMG5_movbdyregpt_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *listv, int ilistv, MMG5_int *lists, int ilists, int improveSurf, int improveVol)
Definition: movpt_3d.c:625
int MMG5_movintpt_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *list, int ilist, int improve)
Definition: movpt_3d.c:58
int MMG3D_curveEndingPts(MMG5_pMesh mesh, MMG5_int *lists, int ilists, const int16_t edgTag, MMG5_int ip0, MMG5_int *ip1, MMG5_int *ip2)
Definition: movpt_3d.c:936
static int MMG3D_movbdycurvept_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *listv, int ilistv, MMG5_int *lists, int ilists, int improve, const int16_t edgTag)
Definition: movpt_3d.c:1255
int MMG3D_movbdycurvept_newPosForSimu(MMG5_pMesh mesh, MMG5_pPoint p0, MMG5_int ip0, MMG5_int ip1, MMG5_int ip2, double ll1old, double ll2old, uint8_t isrid, const double step, double o[3], double no[3], double no2[3], double to[3], const int16_t edgTag)
Definition: movpt_3d.c:1142
int MMG3D_movv_iso(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int k, int ib)
Definition: movpt_3d.c:1875
int MMG5_movintptLES_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, MMG5_int *list, int ilist, int improve)
Definition: movpt_3d.c:185
int MMG3D_movbdyregpt_geom(MMG5_pMesh mesh, MMG5_int *lists, const MMG5_int kel, const MMG5_int ip0, double n[3], double lambda[3], double o[3], double no[3])
Definition: movpt_3d.c:482
int MMG5_movbdyrefpt_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *listv, int ilistv, MMG5_int *lists, int ilists, int improve)
Definition: movpt_3d.c:1438
int MMG5_movbdynompt_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *listv, int ilistv, MMG5_int *lists, int ilists, int improve)
Definition: movpt_3d.c:1462
int8_t parTyp
Definition: libmmgtypes.h:541
double hmax
Definition: libmmgtypes.h:518
MMG5_pPar par
Definition: libmmgtypes.h:517
double hausd
Definition: libmmgtypes.h:518
MMG mesh structure.
Definition: libmmgtypes.h:605
MMG5_Info info
Definition: libmmgtypes.h:651
MMG5_pPoint point
Definition: libmmgtypes.h:641
MMG5_int mark
Definition: libmmgtypes.h:618
MMG5_pxPoint xpoint
Definition: libmmgtypes.h:642
MMG5_int xp
Definition: libmmgtypes.h:620
MMG5_pTetra tetra
Definition: libmmgtypes.h:643
MMG5_pxTetra xtetra
Definition: libmmgtypes.h:644
MMG5_int xpmax
Definition: libmmgtypes.h:612
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
double n[3]
Definition: libmmgtypes.h:272
int16_t tag
Definition: libmmgtypes.h:284
double c[3]
Definition: libmmgtypes.h:271
MMG5_int xp
Definition: libmmgtypes.h:279
MMG5_int ref
Definition: libmmgtypes.h:278
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
MMG5_int ref
Definition: libmmgtypes.h:335
MMG5_int v[3]
Definition: libmmgtypes.h:334
Structure to store surface points of a MMG mesh.
Definition: libmmgtypes.h:294
int8_t nnor
Definition: libmmgtypes.h:297
double n2[3]
Definition: libmmgtypes.h:295
double n1[3]
Definition: libmmgtypes.h:295
Structure to store the surface tetrahedra of a MMG mesh.
Definition: libmmgtypes.h:418
int16_t tag[6]
Definition: libmmgtypes.h:425
int8_t ori
Definition: libmmgtypes.h:427