Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
analys_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"
37#include "libmmg3d_private.h"
38
47 MMG5_pTria ptt;
48 MMG5_int k;
49
50 /* The MG_REQ+MG_NOSURF tag mark the boundary edges that we dont want to touch
51 * but that are not really required (-nosurf option) */
52 for (k=1; k<=mesh->nt; k++) {
53 ptt = &mesh->tria[k];
54
55 if ( mesh->info.nosurf && (!(ptt->tag[0] & MG_REQ)) ) {
56 ptt->tag[0] |= MG_REQ;
57 ptt->tag[0] |= MG_NOSURF;
58 }
59
60 if ( ptt->tag[0] & MG_PARBDY ) {
61 ptt->tag[0] |= MG_NOSURF;
62 ptt->tag[0] |= MG_REQ;
63 }
64
65 if ( mesh->info.nosurf && (!(ptt->tag[1] & MG_REQ)) ) {
66 ptt->tag[1] |= MG_REQ;
67 ptt->tag[1] |= MG_NOSURF;
68 }
69
70 if ( ptt->tag[1] & MG_PARBDY ) {
71 ptt->tag[1] |= MG_NOSURF;
72 ptt->tag[1] |= MG_REQ;
73 }
74
75 if ( mesh->info.nosurf && (!(ptt->tag[2] & MG_REQ)) ) {
76 ptt->tag[2] |= MG_REQ;
77 ptt->tag[2] |= MG_NOSURF;
78 }
79
80 if ( ptt->tag[2] & MG_PARBDY ) {
81 ptt->tag[2] |= MG_NOSURF;
82 ptt->tag[2] |= MG_REQ;
83 }
84 }
85
86 return;
87}
88
89
109 MMG5_pTria pt,pt1;
110 MMG5_int *adja,*adjb,adji1,adji2,*pile,iad,ipil,ip1,ip2,gen;
111 MMG5_int k,kk,iel,jel,nvf,nf,nr,nm,nt,nre,nreq,ncc,ned,ref;
112 int16_t tag;
113 int8_t i,ii,i1,i2,ii1,ii2,voy;
114
115 nvf = nf = ncc = ned = 0;
116
117 MMG5_SAFE_MALLOC(pile,mesh->nt+1,MMG5_int,return 0);
118
119 pile[1] = 1;
120 ipil = 1;
121
122 while ( ipil > 0 ) {
123 ncc++;
124
125 do {
126 k = pile[ipil--];
127 pt = &mesh->tria[k];
128 pt->flag = ncc;
129 if ( !MG_EOK(pt) ) continue;
130
131 adja = &mesh->adjt[3*(k-1)+1];
132 for (i=0; i<3; i++) {
133 if( ((pt->tag[i] & MG_PARBDY) && !(pt->tag[i] & MG_PARBDYBDY)) ||
134 (pt->tag[i] & MG_BDY) ) continue;
135 i1 = MMG5_inxt2[i];
136 i2 = MMG5_iprv2[i];
137 ip1 = pt->v[i1];
138 ip2 = pt->v[i2];
139
140 if ( !mesh->point[ip1].tmp ) mesh->point[ip1].tmp = ++nvf;
141 if ( !mesh->point[ip2].tmp ) mesh->point[ip2].tmp = ++nvf;
142
143 if ( MG_EDG(pt->tag[i]) || pt->tag[i] & MG_REQ ) {
144 tag = mesh->point[ip1].tag;
145 mesh->point[ip1].tag |= pt->tag[i];
146 // Remove the MG_NOSURF tag if the vertex is really required.
147 if ( (tag & MG_REQ) && !(tag & MG_NOSURF) ) {
148 mesh->point[ip1].tag &= ~MG_NOSURF;
149 }
150 tag = mesh->point[ip2].tag;
151 mesh->point[ip2].tag |= pt->tag[i];
152 // Remove the MG_NOSURF tag if the vertex is really required.
153 if ( (tag & MG_REQ) && !(tag & MG_NOSURF) ) {
154 mesh->point[ip2].tag &= ~MG_NOSURF;
155 }
156 }
157
158 /* open boundary */
159 tag = 0;
160 if ( mesh->info.opnbdy ) tag += MG_OPNBDY;
161 if ( !adja[i] ) {
162 /* Mark non-manifold edges and open-boundary ones: note that in open
163 * boundary mode, all non-manifold edges are marked as open (because
164 * at least one of the triangles that shares the edge has no
165 * adjacent). */
166 tag += MG_NOM;
167 pt->tag[i] |= tag;
168 mesh->point[ip1].tag |= tag;
169 mesh->point[ip2].tag |= tag;
170 ned++;
171 continue;
172 }
173
174 kk = adja[i] / 3;
175 ii = adja[i] % 3;
176 if ( kk > k ) ned++;
177
178 /* store adjacent */
179 pt1 = &mesh->tria[kk];
180
181 /* non manifold edge */
182 if ( pt->tag[i] & MG_NOM ) {
183 mesh->point[ip1].tag |= MG_NOM;
184 mesh->point[ip2].tag |= MG_NOM;
185 continue;
186 }
187
188 if ( MMG5_abs(pt1->ref) != MMG5_abs(pt->ref) ) {
189 pt->tag[i] |= MG_REF;
190 pt1->tag[ii] |= MG_REF;
191 mesh->point[ip1].tag |= MG_REF;
192 mesh->point[ip2].tag |= MG_REF;
193 }
194
195 /* store adjacent */
196 if ( !pt1->flag ) {
197 pt1->flag = ncc;
198 pile[++ipil] = kk;
199 }
200
201 /* check orientation */
202 ii1 = MMG5_inxt2[ii];
203 ii2 = MMG5_iprv2[ii];
204 if ( pt1->v[ii1] == ip1 ) {
205 assert ( pt1->base );
206 /* Moebius strip */
207 if ( pt1->base < 0 ) {
208 fprintf(stderr,"\n ## Error: %s: Triangle orientation problem (1):"
209 " Moebius strip?\n",__func__);
210 MMG5_SAFE_FREE(pile);
211 return 0;
212 }
213 /* flip orientation */
214 else {
215 pt1->base = -pt1->base;
216 pt1->v[ii1] = ip2;
217 pt1->v[ii2] = ip1;
218
219 /* update adj */
220 iad = 3*(kk-1)+1;
221 adjb = &mesh->adjt[iad];
222 adji1 = mesh->adjt[iad+ii1];
223 adji2 = mesh->adjt[iad+ii2];
224 adjb[ii1] = adji2;
225 adjb[ii2] = adji1;
226
227 /* modif tag + ref */
228 tag = pt1->tag[ii1];
229 pt1->tag[ii1] = pt1->tag[ii2];
230 pt1->tag[ii2] = tag;
231 ref = pt1->edg[ii1];
232 pt1->edg[ii1] = pt1->edg[ii2];
233 pt1->edg[ii2] = ref;
234
235 /* modif voyeurs */
236 if ( adjb[ii1] ) {
237 iel = adjb[ii1] / 3;
238 voy = adjb[ii1] % 3;
239 mesh->adjt[3*(iel-1)+1+voy] = 3*kk + ii1;
240 }
241 if ( adjb[ii2] ) {
242 iel = adjb[ii2] / 3;
243 voy = adjb[ii2] % 3;
244 mesh->adjt[3*(iel-1)+1+voy] = 3*kk + ii2;
245 }
246 nf++;
247 }
248 }
249 else {
250 /* Mark triangles that have a consistent orientation with their
251 * neighbours */
252 pt1->base = -pt1->base;
253 }
254 }
255 }
256 while ( ipil > 0 );
257
258 /* find next unmarked triangle */
259 ipil = 0;
260 for (kk=1; kk<=mesh->nt; kk++) {
261 pt = &mesh->tria[kk];
262 if ( MG_EOK(pt) && (pt->flag == 0) ) {
263 ipil = 1;
264 pile[ipil] = kk;
265 pt->flag = ncc+1;
266 break;
267 }
268 }
269 }
270
271 /* bilan */
272 nr = nm = nre = nreq = nt = 0;
273 for (k=1; k<=mesh->nt; k++) {
274 pt = &mesh->tria[k];
275 if ( !MG_EOK(pt) ) continue;
276 nt++;
277 adja = &mesh->adjt[3*(k-1)+1];
278 for (i=0; i<3; i++) {
279 if ( ( !MG_EDG(pt->tag[i]) ) && ( !(pt->tag[i] & MG_REQ) ) ) continue;
280
281 jel = adja[i] / 3;
282 if ( !jel || jel > k ) {
283 if ( pt->tag[i] & MG_GEO ) nr++;
284 if ( pt->tag[i] & MG_NOM ) nm++;
285 if ( pt->tag[i] & MG_REF ) nre++;
286 if ( pt->tag[i] & MG_REQ ) nreq++;
287 }
288 }
289 }
290
291 if ( mesh->info.ddebug ) {
292 fprintf(stdout," a- ridges: %" MMG5_PRId " found.\n",nr);
293 fprintf(stdout," a- nm : %" MMG5_PRId " found.\n",nm);
294 fprintf(stdout," a- requir: %" MMG5_PRId " found.\n",nreq);
295 fprintf(stdout," a- connex: %" MMG5_PRId " connected component(s)\n",ncc);
296 fprintf(stdout," a- orient: %" MMG5_PRId " flipped\n",nf);
297 }
298 else if ( abs(mesh->info.imprim) > 3 ) {
299 gen = (2 - nvf + ned - nt) / 2;
300 fprintf(stdout," Connected component: %" MMG5_PRId ", genus: %" MMG5_PRId ", reoriented: %" MMG5_PRId "\n",ncc,gen,nf);
301 fprintf(stdout," Edges: %" MMG5_PRId ", tagged: %" MMG5_PRId ", ridges: %" MMG5_PRId ", required: %" MMG5_PRId ", refs: %" MMG5_PRId "\n",
302 ned,nr+nre+nreq,nr,nreq,nre);
303 }
304
305 MMG5_SAFE_FREE(pile);
306 return 1;
307}
308
311 MMG5_pTria pt,pt1;
312 double n1[3],n2[3],dhd;
313 MMG5_int *adja,k,kk,ne,nr,nrrm;
314 int8_t i,ii,i1,i2;
315 static int8_t warn=0;
316
323 nrrm = 0;
324 for (k=1; k<=mesh->nt; k++) {
325 pt = &mesh->tria[k];
326 if ( !MG_EOK(pt) ) continue;
327
328 /* triangle normal */
329 MMG5_nortri(mesh,pt,n1);
330 adja = &mesh->adjt[3*(k-1)+1];
331 for (i=0; i<3; i++) {
332 if ( ((pt->tag[i] & MG_PARBDY) && !(pt->tag[i] & MG_PARBDYBDY)) ||
333 (pt->tag[i] & MG_BDY) ) continue;
334
335 if ( pt->tag[i] & MG_NOM ) {
336 /* 1. We don't compute ridges along non-manifold edges because: if we
337 * choose to analyze their angle, we have to check the normal deviation
338 * during mesh adaptation (which is not done for now).
339 *
340 * 2. Do not analyze if nm edges are MG_REF ones because we can't
341 * analyze if adjacent tria have same refs due to non-consistency of
342 * adjacency building.
343 */
344 continue;
345 }
346
347 kk = adja[i] / 3;
348 ii = adja[i] % 3;
349
350 if ( kk && k < kk ) {
351 pt1 = &mesh->tria[kk];
352 /* check angle w. neighbor. */
353 MMG5_nortri(mesh,pt1,n2);
354 dhd = n1[0]*n2[0] + n1[1]*n2[1] + n1[2]*n2[2];
355
356 if ( (pt->tag[i] & MG_GEO) || (pt1->tag[ii] & MG_GEO) ) {
357 /* Edge is provided as ridge by the user: check that it isn't at the
358 * interface of 2 triangles belonging to the same plane */
359 if ( fabs(dhd-1.) < MMG5_EPSOK ) {
360 pt->tag[i] &= ~MG_GEO;
361 pt1->tag[ii] &= ~MG_GEO;
362 i1 = MMG5_inxt2[i];
363 i2 = MMG5_inxt2[i1];
364 mesh->point[pt->v[i1]].tag &= ~MG_GEO;
365 mesh->point[pt->v[i2]].tag &= ~MG_GEO;
366 nrrm++;
367 if ( !warn ) {
368 fprintf(stdout,"\n ## Warning: %s: at least one ridge along flat angle.\n"
369 " Ridge tag will be removed from edges but vertices can"
370 " still have tags that interfer with remeshing.\n\n",
371 __func__);
372 warn = 1;
373 }
374 }
375 }
376 }
377 }
378 }
379
382 ne = nr = 0;
383 for (k=1; k<=mesh->nt; k++) {
384 pt = &mesh->tria[k];
385 if ( !MG_EOK(pt) ) continue;
386
387 /* triangle normal */
388 MMG5_nortri(mesh,pt,n1);
389 adja = &mesh->adjt[3*(k-1)+1];
390 for (i=0; i<3; i++) {
391 if ( ((pt->tag[i] & MG_PARBDY) && !(pt->tag[i] & MG_PARBDYBDY)) ||
392 (pt->tag[i] & MG_BDY) ) continue;
393
394 if ( pt->tag[i] & MG_NOM ) {
395 /* 1. We don't compute ridges along non-manifold edges because: if we
396 * choose to analyze their angle, we have to check the normal deviation
397 * during mesh adaptation (which is not done for now).
398 *
399 * 2. Do not analyze if nm edges are MG_REF ones because we can't
400 * analyze if adjacent tria have same refs due to non-consistency of
401 * adjacency building.
402 */
403 continue;
404 }
405
406 kk = adja[i] / 3;
407 ii = adja[i] % 3;
408
409 if ( kk && k < kk ) {
410 pt1 = &mesh->tria[kk];
411 /* reference curve */
412 /* Remark: along non-manifold edges we store only adjacency relationship
413 * between 2 surface parts (other are considered without adja). As we
414 * don't ensure consistency in the choic of the surface we cannot rely
415 * on the current test to detect ref edges. For this reason, all
416 * non-manifold edges have to be marked as reference. */
417 if ( pt1->ref != pt->ref ) {
418 pt->tag[i] |= MG_REF;
419 pt1->tag[ii] |= MG_REF;
420 i1 = MMG5_inxt2[i];
421 i2 = MMG5_inxt2[i1];
422 mesh->point[pt->v[i1]].tag |= MG_REF;
423 mesh->point[pt->v[i2]].tag |= MG_REF;
424 ne++;
425 }
426
427 /* check angle w. neighbor. */
428 MMG5_nortri(mesh,pt1,n2);
429 dhd = n1[0]*n2[0] + n1[1]*n2[1] + n1[2]*n2[2];
430 if ( dhd <= mesh->info.dhd ) {
431 pt->tag[i] |= MG_GEO;
432 pt1->tag[ii] |= MG_GEO;
433 i1 = MMG5_inxt2[i];
434 i2 = MMG5_inxt2[i1];
435 mesh->point[pt->v[i1]].tag |= MG_GEO;
436 mesh->point[pt->v[i2]].tag |= MG_GEO;
437 nr++;
438 }
439 else if ( (pt->tag[i] & MG_GEO) || (pt1->tag[ii] & MG_GEO) ) {
440 /* The MG_GEO tag of a vertex at interface of a "true" ridge and a
441 * "spurious" ridge deleted at step 1 may have been removed by
442 * error */
443 pt->tag[i] |= MG_GEO;
444 pt1->tag[ii] |= MG_GEO;
445 i1 = MMG5_inxt2[i];
446 i2 = MMG5_inxt2[i1];
447 mesh->point[pt->v[i1]].tag |= MG_GEO;
448 mesh->point[pt->v[i2]].tag |= MG_GEO;
449 }
450 }
451 }
452 }
453 if ( abs(mesh->info.imprim) > 3 && nr > 0 ) {
454 fprintf(stdout," %" MMG5_PRId " ridges, %" MMG5_PRId " edges added\n",nr,ne);
455 }
456 if ( abs(mesh->info.imprim) > 3 && nrrm > 0 ) {
457 fprintf(stdout," %" MMG5_PRId " ridges removed\n",nrrm);
458 }
459
460
461 return 1;
462}
463
472 MMG5_pTetra pt;
473 MMG5_pxTetra pxt;
474 MMG5_pPoint ppt;
475 MMG5_int k,lists[MMG3D_LMAX+2];
476 int64_t listv[MMG3D_LMAX+2];
477 int ilists,ilistv;
478 int i0,ier;
479 int8_t i,j;
480 static int8_t mmgWarn = 0;
481
482 for (k=1; k<=mesh->np; k++) {
483 ppt = &mesh->point[k];
484 ppt->s = 0;
485 ppt->flag = mesh->mark;
486 }
487 ++mesh->mark;
488
489 /*count the number of tet around a point*/
490 for (k=1; k<=mesh->ne; k++) {
491 pt = &mesh->tetra[k];
492 if ( !MG_EOK(pt) ) continue;
493
494 for (i=0; i<4; i++) {
495 mesh->point[pt->v[i]].s++;
496 }
497 }
498
499 for (k=1; k<=mesh->ne; k++) {
500 pt = &mesh->tetra[k];
501 if ( !MG_EOK(pt) ) continue;
502
503 /* point j on face i */
504 for (i=0; i<4; i++) {
505 for (j=0; j<3; j++) {
506 if ( pt->xt ) {
507 pxt = &mesh->xtetra[pt->xt];
508 }
509 else pxt = 0;
510
511 i0 = MMG5_idir[i][j];
512 ppt = &mesh->point[pt->v[i0]];
513 if ( !(ppt->tag & MG_BDY) ) continue;
514 if ( ppt->flag == mesh->mark ) continue;
515
516 /* Catch a boundary point by a boundary face */
517 if ( (!pt->xt) || !(MG_BDY & pxt->ftag[i]) ) continue;
518 if( ppt->tag & MG_NOM ){
519 if ( mesh->adja[4*(k-1)+1+i] ) continue;
520 ier=MMG5_boulesurfvolp(mesh,k,i0,i,listv,&ilistv,lists,&ilists,1);
521 } else {
522 ier=MMG5_boulesurfvolp(mesh,k,i0,i,listv,&ilistv,lists,&ilists,0);
523 }
524 if ( ier != 1 && !mmgWarn ) {
525 mmgWarn = 1;
526 printf(" ## Warning: %s: unable to check that we don't have"
527 " non-connected domains.\n",__func__);
528 }
529
530 if(ilistv != ppt->s) {
531 if(!(ppt->tag & MG_REQ) ) {
532 ppt->tag |= MG_REQ;
533 ppt->tag |= MG_CRN;
534 }
535 }
536 ppt->flag = mesh->mark;
537 }
538 }
539 }
540 return 1;
541}
542
545 MMG5_pTria pt;
546 MMG5_pPoint ppt,p1,p2;
547 double ux,uy,uz,vx,vy,vz,dd;
548 MMG5_int list[MMG3D_LMAX+2],listref[MMG3D_LMAX+2],k,nc,nre;
549 int xp,nr,ns;
550 int8_t i;
551
552 nre = nc = 0;
553 for (k=1; k<=mesh->nt; k++) {
554 pt = &mesh->tria[k];
555 if ( !MG_EOK(pt) ) continue;
556
557 for (i=0; i<3; i++) {
558 ppt = &mesh->point[pt->v[i]];
559 if ( !MG_VOK(ppt) || ( ppt->tag & MG_CRN ) || ( ppt->tag & MG_NOM ) ||
560 ( ppt->tag & MG_PARBDY ) ) continue;
561 else if ( MG_EDG(ppt->tag) ) {
562 /* Store the number of ridges passing through the point (xp) and the
563 * number of ref edges (nr) */
564 ns = MMG5_bouler(mesh,mesh->adjt,k,i,list,listref,&xp,&nr,MMG3D_LMAX);
565
566 if ( !ns ) continue;
567 if ( (xp+nr) > 2 ) {
568 ppt->tag |= MG_CRN + MG_REQ;
569 ppt->tag &= ~MG_NOSURF;
570 nre++;
571 nc++;
572 }
573 else if ( (xp == 1) && (nr == 1) ) {
574 ppt->tag |= MG_REQ;
575 ppt->tag &= ~MG_NOSURF;
576 nre++;
577 }
578 else if ( xp == 1 && !nr ){
579 ppt->tag |= MG_CRN + MG_REQ;
580 ppt->tag &= ~MG_NOSURF;
581 nre++;
582 nc++;
583 }
584 else if ( nr == 1 && !xp ){
585 ppt->tag |= MG_CRN + MG_REQ;
586 ppt->tag &= ~MG_NOSURF;
587 nre++;
588 nc++;
589 }
590 /* check ridge angle */
591 else {
592 p1 = &mesh->point[list[1]];
593 p2 = &mesh->point[list[2]];
594 ux = p1->c[0] - ppt->c[0];
595 uy = p1->c[1] - ppt->c[1];
596 uz = p1->c[2] - ppt->c[2];
597 vx = p2->c[0] - ppt->c[0];
598 vy = p2->c[1] - ppt->c[1];
599 vz = p2->c[2] - ppt->c[2];
600 dd = (ux*ux + uy*uy + uz*uz) * (vx*vx + vy*vy + vz*vz);
601 if ( fabs(dd) > MMG5_EPSD ) {
602 dd = (ux*vx + uy*vy + uz*vz) / sqrt(dd);
603 if ( dd > -mesh->info.dhd ) {
604 ppt->tag |= MG_CRN;
605 nc++;
606 }
607 }
608 }
609 }
610 }
611 }
612
613 if ( abs(mesh->info.imprim) > 3 && nre > 0 )
614 fprintf(stdout," %" MMG5_PRId " corners, %" MMG5_PRId " singular points detected\n",nc,nre);
615 return 1;
616}
617
639 MMG5_pTria pt;
640 MMG5_pPoint ppt;
641 MMG5_xPoint *pxp;
642 double n[3],dd;
643 MMG5_int *adja,k,kk,ng,nn,nt,nf,nnr;
644 int8_t i,ii,i1;
645
647 if ( mesh->xpoint ) {
648 if ( abs(mesh->info.imprim) > 3 || mesh->info.ddebug ) {
649 fprintf(stdout," ## Warning: %s: no research of boundary points"
650 " and normals of mesh. mesh->xpoint must be freed to enforce"
651 " analysis.\n",__func__);
652 }
653 return 1;
654 }
655
657 ++mesh->base;
658 mesh->xp = 0;
659 nnr = 0;
660 for (k=1; k<=mesh->nt; k++) {
661 pt = &mesh->tria[k];
662 if ( !MG_EOK(pt) ) continue;
663
664 for (i=0; i<3; i++) {
665 ppt = &mesh->point[pt->v[i]];
666 if ( ppt->flag == mesh->base ) continue;
667 else {
668 ++mesh->xp;
669 ppt->flag = mesh->base;
670 if ( mesh->nc1 ) {
671 if ( ppt->n[0]*ppt->n[0] + ppt->n[1]*ppt->n[1] + ppt->n[2]*ppt->n[2] > 0 ) {
676 if ( ppt->tag & MG_PARBDY || ppt->tag & MG_CRN || MG_EDG_OR_NOM(ppt->tag) ) {
677 ++nnr;
678 continue;
679 }
680 ppt->xp = -1;
681 }
682 }
683 }
684 }
685 }
686
688 mesh->xpmax = MG_MAX( (long long)(1.5*mesh->xp),mesh->npmax);
689
690 MMG5_ADD_MEM(mesh,(mesh->xpmax+1)*sizeof(MMG5_xPoint),"boundary points",return 0);
692
694 nn = ng = nt = nf = 0;
695 mesh->xp = 0;
696 ++mesh->base;
697 for (k=1; k<=mesh->nt; k++) {
698 pt = &mesh->tria[k];
699 if ( !MG_EOK(pt) ) continue;
700
701 adja = &mesh->adjt[3*(k-1)+1];
702 for (i=0; i<3; i++) {
703 ppt = &mesh->point[pt->v[i]];
704
705 if ( ppt->tag & MG_PARBDY || ppt->tag & MG_CRN || ppt->tag & MG_NOM || ppt->flag == mesh->base )
706 {
707 continue;
708 }
709
711 if ( !MG_EDG(ppt->tag) ) {
712
713 if ( (!mesh->nc1) ||
714 ppt->n[0]*ppt->n[0]+ppt->n[1]*ppt->n[1]+ppt->n[2]*ppt->n[2]<=MMG5_EPSD2 ) {
715 if ( !MMG5_boulen(mesh,mesh->adjt,k,i,ppt->n) ) {
716 ++nf;
717 continue;
718 }
719 else ++nn;
720 }
721
722 ++mesh->xp;
723 if(mesh->xp > mesh->xpmax){
725 "larger xpoint table",
726 mesh->xp--;return 0;);
727 }
728 ppt->xp = mesh->xp;
729 pxp = &mesh->xpoint[ppt->xp];
730 memcpy(pxp->n1,ppt->n,3*sizeof(double));
731 ppt->n[0] = ppt->n[1] = ppt->n[2] = 0.;
732 ppt->flag = mesh->base;
733
734 }
735
737 i1 = MMG5_inxt2[i];
738 if ( !MG_EDG(pt->tag[i1]) ) {
739 continue;
740 }
741
742 /* As we skip non-manifold point, the edge should be manifold */
743 assert ( (!(MG_NOM & pt->tag[i1])) && "Unexpected non-manifold edge" );
744 if ( !MMG5_boulen(mesh,mesh->adjt,k,i,n) ) {
745 ++nf;
746 continue;
747 }
748
749 ++mesh->xp;
750 if(mesh->xp > mesh->xpmax){
752 "larger xpoint table",
753 mesh->xp--;return 0;);
754 }
755 ppt->xp = mesh->xp;
756 pxp = &mesh->xpoint[ppt->xp];
757 memcpy(pxp->n1,n,3*sizeof(double));
758
759 if ( (pt->tag[i1] & MG_GEO) && adja[i1] > 0 ) {
760 kk = adja[i1] / 3;
761 ii = adja[i1] % 3;
762 ii = MMG5_inxt2[ii];
763 if ( !MMG5_boulen(mesh,mesh->adjt,kk,ii,n) ) {
764 ++nf;
765 continue;
766 }
767 memcpy(pxp->n2,n,3*sizeof(double));
768
770 ppt->n[0] = pxp->n1[1]*pxp->n2[2] - pxp->n1[2]*pxp->n2[1];
771 ppt->n[1] = pxp->n1[2]*pxp->n2[0] - pxp->n1[0]*pxp->n2[2];
772 ppt->n[2] = pxp->n1[0]*pxp->n2[1] - pxp->n1[1]*pxp->n2[0];
773 dd = ppt->n[0]*ppt->n[0] + ppt->n[1]*ppt->n[1] + ppt->n[2]*ppt->n[2];
774 if ( dd > MMG5_EPSD2 ) {
775 dd = 1.0 / sqrt(dd);
776 ppt->n[0] *= dd;
777 ppt->n[1] *= dd;
778 ppt->n[2] *= dd;
779 }
780 ppt->flag = mesh->base;
781 ++nt;
782 continue;
783 }
784
785 /* compute tgte */
786 ppt->flag = mesh->base;
787 ++nt;
788 if ( !MMG5_boulec(mesh,mesh->adjt,k,i,ppt->n) ) {
789 ++nf;
790 continue;
791 }
792 dd = pxp->n1[0]*ppt->n[0] + pxp->n1[1]*ppt->n[1] + pxp->n1[2]*ppt->n[2];
793 ppt->n[0] -= dd*pxp->n1[0];
794 ppt->n[1] -= dd*pxp->n1[1];
795 ppt->n[2] -= dd*pxp->n1[2];
796 dd = ppt->n[0]*ppt->n[0] + ppt->n[1]*ppt->n[1] + ppt->n[2]*ppt->n[2];
797 if ( dd > MMG5_EPSD2 ) {
798 dd = 1.0 / sqrt(dd);
799 ppt->n[0] *= dd;
800 ppt->n[1] *= dd;
801 ppt->n[2] *= dd;
802 }
803 }
804 }
805 mesh->nc1 = 0;
806
807 if ( abs(mesh->info.imprim) > 3 && nn+nt > 0 ) {
808 if ( nnr )
809 fprintf(stdout," %" MMG5_PRId " input normals ignored\n",nnr);
810 fprintf(stdout," %" MMG5_PRId " normals, %" MMG5_PRId " tangents updated (%" MMG5_PRId " failed)\n",nn,nt,nf);
811 }
812 return 1;
813}
814
828static inline int MMG3D_dichotomytria(MMG5_pMesh mesh, MMG5_pTria pt, MMG5_int k, double *c, double *n) {
829
830 MMG5_pPoint ppt;
831 double to,tp,t,nnew[3],p[3],o[3],result;
832 int it,maxit,pos,ier;
833
834 it = 0;
835 maxit = 5;
836 to = 0.0;
837 tp = 1.0;
838 t = 0.5;
839 pos = 0;
840
841 ppt = &mesh->point[k];
842
843 /* initial coordinates of point before regularization */
844 o[0] = ppt->c[0];
845 o[1] = ppt->c[1];
846 o[2] = ppt->c[2];
847
848 /* initial coordinates of new point */
849 p[0] = c[0];
850 p[1] = c[1];
851 p[2] = c[2];
852
853 do {
854 mesh->point[0].c[0] = o[0] + t*(p[0] - o[0]);
855 mesh->point[0].c[1] = o[1] + t*(p[1] - o[1]);
856 mesh->point[0].c[2] = o[2] + t*(p[2] - o[2]);
857
858 MMG5_nortri(mesh, pt, nnew);
859 MMG5_dotprod(3,n,nnew,&result);
860
861 if ( result <= 0.0 ) {
862 tp = t;
863 }
864 else {
865 c[0] = mesh->point[0].c[0];
866 c[1] = mesh->point[0].c[1];
867 c[2] = mesh->point[0].c[2];
868 to = t;
869 pos = 1;
870 }
871
872 t = 0.5*(to + tp);
873 }
874 while ( ++it < maxit );
875
876 if ( pos ) {
877 return 1;
878 }
879 else
880 return 0;
881}
882
895static inline int MMG3D_dichotomytetra(MMG5_pMesh mesh, MMG5_int *v, MMG5_int k, double *c) {
896
897 MMG5_pPoint ppt;
898 double p[3],o[3],vol,to,tp,t;
899 int it,maxit,pos;
900
901 it = 0;
902 maxit = 5;
903 to = 0.0;
904 tp = 1.0;
905 t = 0.5;
906 pos = 0;
907
908 ppt = &mesh->point[k];
909
910 /* initial coordinates of point before regularization */
911 o[0] = ppt->c[0];
912 o[1] = ppt->c[1];
913 o[2] = ppt->c[2];
914
915 /* initial coordinates of new point */
916 p[0] = c[0];
917 p[1] = c[1];
918 p[2] = c[2];
919
920 do {
921 mesh->point[0].c[0] = o[0] + t*(p[0] - o[0]);
922 mesh->point[0].c[1] = o[1] + t*(p[1] - o[1]);
923 mesh->point[0].c[2] = o[2] + t*(p[2] - o[2]);
924
925 vol = MMG5_orvol(mesh->point,v);
926
927 if ( vol <= 0.0 ) {
928 tp = t;
929 }
930 else {
931 c[0] = mesh->point[0].c[0];
932 c[1] = mesh->point[0].c[1];
933 c[2] = mesh->point[0].c[2];
934 to = t;
935 pos = 1;
936 }
937
938 t = 0.5*(to + tp);
939 }
940 while ( ++it < maxit );
941
942 if ( pos ) {
943 return 1;
944 }
945 else
946 return 0;
947}
948
957 MMG5_pTria pt;
958 MMG5_pTetra ptet;
959 MMG5_pPoint ppt,p0;
960 MMG5_Tria tnew;
961 double *tabl,c[3],n[3],nnew[3],*cptr,lm1,lm2,cx,cy,cz,res0,res,result;
962 int i,ii,it,nit,ilist,noupdate;
963 MMG5_int k,kt,nn,iel,list[MMG5_LMAX],tlist[MMG5_LMAX],*adja,iad,v[4];
964 int64_t tetlist[MMG5_LMAX];
965
966 /* assign seed to vertex */
967 for (k=1; k<=mesh->nt; k++) {
968 pt = &mesh->tria[k];
969 if ( !MG_EOK(pt) ) continue;
970 for (i=0; i<3; i++) {
971 ppt = &mesh->point[pt->v[i]];
972 if ( !ppt->s ) ppt->s = k;
973 }
974 }
975
976 for (k=1; k<=mesh->np; k++) {
977 ppt = &mesh->point[k];
978 ppt->flag = 0;
979 }
980
981 for (k=1; k<=mesh->ne; k++) {
982 ptet = &mesh->tetra[k];
983 if ( !MG_EOK(ptet) ) continue;
984 for (i=0; i<4; i++) {
985 ppt = &mesh->point[ptet->v[i]];
986 if ( !ppt->flag ) ppt->flag = k;
987 }
988 }
989
990 /* allocate memory for coordinates */
991 MMG5_SAFE_CALLOC(tabl,3*mesh->np+1,double,return 0);
992
993 /* Pointer toward the suitable adjacency array for Mmgs and Mmg3d */
994 adja = mesh->adjt;
995
996 it = 0;
997 nit = 10;
998 res0 = 0.0;
999 lm1 = 0.4;
1000 lm2 = 0.399;
1001 while ( it++ < nit ) {
1002 /* step 1: laplacian */
1003 for (k=1; k<=mesh->np; k++) {
1004 ppt = &mesh->point[k];
1005
1006 iad = 3*(k-1)+1;
1007 tabl[iad+0] = ppt->c[0];
1008 tabl[iad+1] = ppt->c[1];
1009 tabl[iad+2] = ppt->c[2];
1010 if ( !MG_VOK(ppt) ) continue;
1011 if ( ppt->tag & MG_CRN || ppt->tag & MG_NOM || MG_EDG(ppt->tag) ) continue;
1012
1013 iel = ppt->s;
1014 if ( !iel ) continue; // Mmg3d
1015
1016 pt = &mesh->tria[iel];
1017 i = 0;
1018 if ( pt->v[1] == k ) i = 1;
1019 else if ( pt->v[2] == k ) i = 2;
1020
1021 ilist = MMG5_boulep(mesh,iel,i,adja,list,tlist);
1022
1023 /* average coordinates */
1024 cx = cy = cz = 0.0;
1025 for (i=1; i<=ilist; i++) {
1026 p0 = &mesh->point[list[i]];
1027
1028 cptr = p0->c;
1029 cx += cptr[0];
1030 cy += cptr[1];
1031 cz += cptr[2];
1032 }
1033 cx /= ilist;
1034 cy /= ilist;
1035 cz /= ilist;
1036
1037 /* Laplacian */
1038 cptr = ppt->c;
1039 tabl[iad+0] = cptr[0] + lm1 * (cx - cptr[0]);
1040 tabl[iad+1] = cptr[1] + lm1 * (cy - cptr[1]);
1041 tabl[iad+2] = cptr[2] + lm1 * (cz - cptr[2]);
1042 }
1043
1044 /* step 2: anti-laplacian */
1045 res = 0;
1046 nn = 0;
1047 for (k=1; k<=mesh->np; k++) {
1048 ppt = &mesh->point[k];
1049
1050 if ( !MG_VOK(ppt) ) continue;
1051 if ( ppt->tag & MG_CRN || ppt->tag & MG_NOM || MG_EDG(ppt->tag) ) continue;
1052
1053 iel = ppt->s;
1054 if ( !iel ) continue; // Mmg3d
1055
1056 pt = &mesh->tria[iel];
1057 i = 0;
1058 if ( pt->v[1] == k ) i = 1;
1059 else if ( pt->v[2] == k ) i = 2;
1060
1061 ilist = MMG5_boulep(mesh,iel,i,adja,list,tlist);
1062
1063 /* average normal */
1064 cx = cy = cz = 0.0;
1065 for (i=1; i<=ilist; i++) {
1066 iad = 3*(list[i]-1) + 1;
1067 cx += tabl[iad+0];
1068 cy += tabl[iad+1];
1069 cz += tabl[iad+2];
1070 }
1071 cx /= ilist;
1072 cy /= ilist;
1073 cz /= ilist;
1074
1075 /* antiLaplacian */
1076 iad = 3*(k-1)+1;
1077 c[0] = tabl[iad+0] - lm2 * (cx - tabl[iad+0]);
1078 c[1] = tabl[iad+1] - lm2 * (cy - tabl[iad+1]);
1079 c[2] = tabl[iad+2] - lm2 * (cz - tabl[iad+2]);
1080
1081 cptr = ppt->c;
1082
1083 mesh->point[0].c[0] = c[0];
1084 mesh->point[0].c[1] = c[1];
1085 mesh->point[0].c[2] = c[2];
1086
1087 /* check for negative areas */
1088 noupdate = 0;
1089 for (kt = 0 ; kt<ilist ; kt++) {
1090 pt = &mesh->tria[tlist[kt]];
1091
1092 if ( !MG_EOK(pt) ) continue;
1093
1094 MMG5_nortri(mesh, pt, n);
1095
1096 for (i=0;i<3;i++) {
1097 tnew.v[i] = pt->v[i];
1098 }
1099
1100 i = 0;
1101 if ( pt->v[1] == k ) i = 1;
1102 if ( pt->v[2] == k ) i = 2;
1103
1104 tnew.v[i] = 0;
1105
1106 MMG5_nortri(mesh, &tnew, nnew);
1107 MMG5_dotprod(3,n,nnew,&result);
1108 if ( result < 0.0 ) {
1109 if (!MMG3D_dichotomytria(mesh,&tnew,k,c,n))
1110 noupdate = 1;
1111 continue;
1112 }
1113 }
1114
1115 /* check for negative volumes */
1116 if ( !noupdate ) {
1117 iel = ppt->flag;
1118 ptet = &mesh->tetra[iel];
1119
1120 if ( !MG_EOK(ptet) ) continue;
1121
1122 i = 0;
1123 if ( ptet->v[1] == k ) i = 1;
1124 else if ( ptet->v[2] == k ) i = 2;
1125 else if ( ptet->v[3] == k ) i = 3;
1126 ilist = MMG5_boulevolp(mesh, iel, i, tetlist);
1127
1128 for( kt=0 ; kt<ilist ; kt++ ) {
1129 iel = tetlist[kt] / 4;
1130 i = tetlist[kt] % 4;
1131 ptet = &mesh->tetra[iel];
1132
1133 for ( ii=0 ; ii<4 ; ii++ ) {
1134 v[ii] = ptet->v[ii];
1135 }
1136 v[i] = 0;
1137 result = MMG5_orvol(mesh->point,v);
1138
1139 if ( result <= 0.0 ) {
1140 if (!MMG3D_dichotomytetra(mesh,v,k,c))
1141 noupdate = 1;
1142 continue;
1143 }
1144 }
1145 }
1146 if ( !noupdate ) {
1147 res += (cptr[0]-c[0])*(cptr[0]-c[0]) + (cptr[1]-c[1])*(cptr[1]-c[1]) + (cptr[2]-c[2])*(cptr[2]-c[2]);
1148 cptr[0] = c[0];
1149 cptr[1] = c[1];
1150 cptr[2] = c[2];
1151 nn++;
1152 }
1153 }
1154 if ( it == 1 ) res0 = res;
1155 if ( res0 > MMG5_EPSD ) res = res / res0;
1156 if ( mesh->info.imprim < -1 || mesh->info.ddebug ) {
1157 fprintf(stdout," iter %5d res %.3E\r",it,res);
1158 fflush(stdout);
1159 }
1160 if ( it > 1 && res < MMG5_EPS ) break;
1161 }
1162 /* reset the ppt->s and ppt->flag tags */
1163 for (k=1; k<=mesh->np; ++k) {
1164 mesh->point[k].s = 0;
1165 mesh->point[k].flag = 0;
1166 }
1167
1168 if ( mesh->info.imprim < -1 || mesh->info.ddebug ) fprintf(stdout,"\n");
1169
1170 if ( abs(mesh->info.imprim) > 4 )
1171 fprintf(stdout," %" MMG5_PRId " coordinates regularized: %.3e\n",nn,res);
1172
1173 MMG5_SAFE_FREE(tabl);
1174 return 1;
1175}
1176
1177
1199 MMG5_pTetra pt;
1200 MMG5_pPoint p0;
1201 MMG5_pxPoint pxp;
1202 MMG5_int k;
1203 MMG5_int *adja,base;
1204 double n[3],t[3];
1205 int8_t i,j,ip,ier;
1206
1207 base = ++mesh->base;
1208 for (k=1; k<=mesh->ne; k++) {
1209 pt = &mesh->tetra[k];
1210 if( !MG_EOK(pt) ) continue;
1211 adja = &mesh->adja[4*(k-1)+1];
1212 for (i=0; i<4; i++) {
1213 if ( adja[i] ) continue;
1214 for (j=0; j<3; j++) {
1215 ip = MMG5_idir[i][j];
1216 p0 = &mesh->point[pt->v[ip]];
1217 if ( p0->flag == base ) continue;
1218 else if ( (!(p0->tag & MG_NOM)) || (p0->tag & MG_PARBDY) ) {
1219 continue;
1220 }
1221
1222 p0->flag = base;
1223 ier = MMG5_boulenm(mesh,k,ip,i,n,t);
1224
1225 if ( ier < 0 )
1226 return 0;
1227 else if ( !ier ) {
1228 p0->tag |= MG_REQ;
1229 p0->tag &= ~MG_NOSURF;
1230 }
1231 else {
1232 if ( !p0->xp ) {
1233 ++mesh->xp;
1234 if(mesh->xp > mesh->xpmax){
1236 "larger xpoint table",
1237 mesh->xp--;
1238 fprintf(stderr," Exit program.\n");return 0;);
1239 }
1240 p0->xp = mesh->xp;
1241 }
1242 pxp = &mesh->xpoint[p0->xp];
1243 memcpy(pxp->n1,n,3*sizeof(double));
1244 memcpy(p0->n,t,3*sizeof(double));
1245 }
1246 }
1247 }
1248 }
1249 /* Deal with the non-manifold points that do not belong to a surface
1250 * tetra (a tetra that has a face without adjacent)*/
1251 for (k=1; k<=mesh->ne; k++) {
1252 pt = &mesh->tetra[k];
1253 if( !MG_EOK(pt) ) continue;
1254
1255 for (i=0; i<4; i++) {
1256 p0 = &mesh->point[pt->v[i]];
1257 if ( p0->tag & MG_REQ || !(p0->tag & MG_NOM) || p0->xp || (p0->tag & MG_PARBDY) ) {
1258 /* CRN points are skipped because we skip REQ points.
1259 * For connex meshes it is possible
1260 * to compute the tangent at required points but it is not if the mesh
1261 * contains some edge or point-connections: all points along such a feature line are marked
1262 * as required by the previous loop because \ref MMB5_boulnm fails, thus we pass here
1263 * even if the line belongs to a surface. Then \ref MMG5_boulenmInt fails too and
1264 * we end up without normals and tangents. For sake of simplicity and because normals
1265 * and tangents at required points are not used for remahsing, we choose to skip all
1266 * internal REQ points. */
1267 continue;
1268 }
1269 ier = MMG5_boulenmInt(mesh,k,i,t);
1270 if ( ier ) {
1271 ++mesh->xp;
1272 if(mesh->xp > mesh->xpmax){
1274 "larger xpoint table",
1275 mesh->xp--;
1276 fprintf(stderr," Exit program.\n");return 0;);
1277 }
1278 p0->xp = mesh->xp;
1279 pxp = &mesh->xpoint[p0->xp];
1280 memcpy(p0->n,t,3*sizeof(double));
1281 pxp->nnor = 1;
1282 }
1283 else {
1284 p0->tag |= MG_REQ;
1285 p0->tag &= ~MG_NOSURF;
1286 }
1287 }
1288 }
1289
1290 /*for (k=1; k<=mesh->np; k++) {
1291 p0 = &mesh->point[k];
1292 if ( !(p0->tag & MG_NOM) || p0->xp ) continue;
1293 p0->tag |= MG_REQ;
1294 p0->tag &= ~MG_NOSURF;
1295 }*/
1296
1297 return 1;
1298}
1299
1325 MMG5_Hash hash;
1326 int ier;
1327
1329 if ( abs(mesh->info.imprim) > 3 )
1330 fprintf(stdout,"\n ** SURFACE ANALYSIS\n");
1331
1332 /* create tetra adjacency */
1333 if ( !MMG3D_hashTetra(mesh,1) ) {
1334 fprintf(stderr,"\n ## Hashing problem (1). Exit program.\n");
1335 return 0;
1336 }
1337
1338 if ( mesh->info.iso && mesh->info.opnbdy ) {
1340 if ( !ier ) {
1341 fprintf(stderr,"\n ## Problem when updating the xtetra data after ls discretization."
1342 " Exit program.\n");
1343 return 0;
1344 }
1345 }
1346
1347 /* create prism adjacency */
1348 if ( !MMG3D_hashPrism(mesh) ) {
1349 fprintf(stderr,"\n ## Prism hashing problem. Exit program.\n");
1350 return 0;
1351 }
1352
1353 /* compatibility triangle orientation w/r tetras */
1354 if ( !MMG5_bdryPerm(mesh) ) {
1355 fprintf(stderr,"\n ## Boundary orientation problem. Exit program.\n");
1356 return 0;
1357 }
1358
1359 /* identify surface mesh */
1360 if ( !MMG5_chkBdryTria(mesh) ) {
1361 fprintf(stderr,"\n ## Boundary problem. Exit program.\n");
1362 return 0;
1363 }
1364
1367
1368 /* Set surface triangles to required in nosurf mode or for parallel boundaries */
1370
1371 /* create surface adjacency */
1372 memset ( &hash, 0x0, sizeof(MMG5_Hash));
1373 if ( !MMG3D_hashTria(mesh,&hash) ) {
1374 MMG5_DEL_MEM(mesh,hash.item);
1375 fprintf(stderr,"\n ## Hashing problem (2). Exit program.\n");
1376 return 0;
1377 }
1378
1379 /* build hash table for geometric edges */
1380 if ( !MMG5_hGeom(mesh) ) {
1381 fprintf(stderr,"\n ## Hashing problem (0). Exit program.\n");
1382 MMG5_DEL_MEM(mesh,hash.item);
1384 return 0;
1385 }
1386
1388 if ( abs(mesh->info.imprim) > 5 || mesh->info.ddebug )
1389 fprintf(stdout," ** SETTING TOPOLOGY\n");
1390
1391 /* identify connexity, flip orientation of faces if needed and transfer
1392 * triangle edge tags toward vertices. */
1393 if ( !MMG5_setadj(mesh) ) {
1394 fprintf(stderr,"\n ## Topology problem. Exit program.\n");
1395 MMG5_DEL_MEM(mesh,hash.item);
1396 return 0;
1397 }
1398
1399 /* check for ridges */
1400 if ( mesh->info.dhd > MMG5_ANGLIM && !MMG5_setdhd(mesh) ) {
1401 fprintf(stderr,"\n ## Geometry problem. Exit program.\n");
1402 MMG5_DEL_MEM(mesh,hash.item);
1403 return 0;
1404 }
1405
1406 /* identify singularities */
1407 if ( !MMG5_singul(mesh) ) {
1408 fprintf(stderr,"\n ## MMG5_Singularity problem. Exit program.\n");
1409 MMG5_DEL_MEM(mesh,hash.item);
1410 return 0;
1411 }
1412
1413 if ( abs(mesh->info.imprim) > 3 || mesh->info.ddebug )
1414 fprintf(stdout," ** DEFINING GEOMETRY\n");
1415
1416 /* regularize vertices coordinates*/
1417 if ( mesh->info.xreg && !MMG3D_regver(mesh) ) {
1418 fprintf(stderr,"\n ## Coordinates regularization problem. Exit program.\n");
1419 return 0;
1420 }
1421
1422 /* define (and regularize) normals */
1423 if ( !MMG5_norver(mesh) ) {
1424 fprintf(stderr,"\n ## Normal problem. Exit program.\n");
1425 MMG5_DEL_MEM(mesh,hash.item);
1426 return 0;
1427 }
1428 if ( mesh->info.nreg && !MMG5_regnor(mesh) ) {
1429 fprintf(stderr,"\n ## Normal regularization problem. Exit program.\n");
1430 return 0;
1431 }
1432
1433 /* set bdry entities to tetra and fill the orientation field */
1434 if ( !MMG5_bdrySet(mesh) ) {
1435 fprintf(stderr,"\n ## Boundary problem. Exit program.\n");
1436 MMG5_DEL_MEM(mesh,hash.item);
1438 return 0;
1439 }
1440
1441 /* set non-manifold edges sharing non-intersecting multidomains as required */
1442 if ( abs(mesh->info.imprim) > 5 || mesh->info.ddebug )
1443 fprintf(stdout," ** UPDATING TOPOLOGY AT NON-MANIFOLD POINTS\n");
1444
1445 if ( !MMG5_setNmTag(mesh,&hash) ) {
1446 fprintf(stderr,"\n ## Non-manifold topology problem. Exit program.\n");
1447 MMG5_DEL_MEM(mesh,hash.item);
1449 return 0;
1450 }
1451
1452 /* check subdomains connected by a vertex and mark these vertex as corner and
1453 required */
1455
1456 /* build hash table for geometric edges */
1457 if ( !mesh->na && !MMG5_hGeom(mesh) ) {
1458 fprintf(stderr,"\n ## Hashing problem (0). Exit program.\n");
1461 return 0;
1462 }
1463
1464 /* Update edges tags and references for xtetras */
1465 if ( !MMG5_bdryUpdate(mesh) ) {
1466 fprintf(stderr,"\n ## Boundary problem. Exit program.\n");
1468 return 0;
1469 }
1470
1471 /* define geometry for non manifold points */
1472 if ( !MMG3D_nmgeom(mesh) ) return 0;
1473
1474#ifdef USE_POINTMAP
1475 /* Initialize source point with input index */
1476 MMG5_int ip;
1477 for( ip = 1; ip <= mesh->np; ip++ )
1478 mesh->point[ip].src = ip;
1479#endif
1480
1481 /* release memory */
1485 mesh->nt = 0;
1486
1488
1489 return 1;
1490}
int ier
MMG5_pMesh * mesh
int MMG5_regnor(MMG5_pMesh mesh)
Definition: analys.c:46
void MMG3D_set_reqBoundaries(MMG5_pMesh mesh)
Definition: analys_3d.c:46
int MMG3D_analys(MMG5_pMesh mesh)
Definition: analys_3d.c:1324
int MMG5_setdhd(MMG5_pMesh mesh)
Definition: analys_3d.c:310
int MMG5_norver(MMG5_pMesh mesh)
Definition: analys_3d.c:638
static int MMG3D_dichotomytetra(MMG5_pMesh mesh, MMG5_int *v, MMG5_int k, double *c)
Definition: analys_3d.c:895
int MMG3D_regver(MMG5_pMesh mesh)
Definition: analys_3d.c:956
static int MMG3D_dichotomytria(MMG5_pMesh mesh, MMG5_pTria pt, MMG5_int k, double *c, double *n)
Definition: analys_3d.c:828
int MMG5_setadj(MMG5_pMesh mesh)
Definition: analys_3d.c:108
int MMG5_chkVertexConnectedDomains(MMG5_pMesh mesh)
Definition: analys_3d.c:471
int MMG5_singul(MMG5_pMesh mesh)
Definition: analys_3d.c:544
int MMG3D_nmgeom(MMG5_pMesh mesh)
Definition: analys_3d.c:1198
int MMG5_boulec(MMG5_pMesh mesh, MMG5_int *adjt, MMG5_int start, int ip, double *tt)
Definition: boulep.c:199
MMG5_Info info
int MMG5_boulep(MMG5_pMesh mesh, MMG5_int start, int ip, MMG5_int *adja, MMG5_int *list, MMG5_int *tlist)
Definition: boulep.c:52
int MMG5_bouler(MMG5_pMesh mesh, MMG5_int *adjt, MMG5_int start, int ip, MMG5_int *list, MMG5_int *listref, int *ng, int *nr, int lmax)
Definition: boulep.c:287
int MMG5_boulen(MMG5_pMesh mesh, MMG5_int *adjt, MMG5_int start, int ip, double *nn)
Definition: boulep.c:119
int MMG5_boulevolp(MMG5_pMesh mesh, MMG5_int start, int ip, int64_t *list)
Definition: boulep_3d.c:54
int MMG5_boulesurfvolp(MMG5_pMesh mesh, MMG5_int start, int ip, int iface, int64_t *listv, int *ilistv, MMG5_int *lists, int *ilists, int isnm)
Definition: boulep_3d.c:607
int MMG5_boulenm(MMG5_pMesh mesh, MMG5_int start, int ip, int iface, double n[3], double t[3])
Definition: boulep_3d.c:196
int MMG5_boulenmInt(MMG5_pMesh mesh, MMG5_int start, int ip, double t[3])
Definition: boulep_3d.c:342
#define MMG5_EPSD
#define MMG5_EPS
int MMG3D_hashTetra(MMG5_pMesh mesh, int pack)
Definition: hash_3d.c:122
int MMG5_hGeom(MMG5_pMesh mesh)
Definition: hash_3d.c:1022
int MMG3D_hashPrism(MMG5_pMesh mesh)
Definition: hash_3d.c:238
int MMG5_chkBdryTria(MMG5_pMesh mesh)
Definition: hash_3d.c:1439
int MMG5_bdrySet(MMG5_pMesh mesh)
Definition: hash_3d.c:1744
int MMG5_bdryUpdate(MMG5_pMesh mesh)
Definition: hash_3d.c:2044
int MMG5_setNmTag(MMG5_pMesh mesh, MMG5_Hash *hash)
Definition: hash_3d.c:727
int MMG3D_hashTria(MMG5_pMesh mesh, MMG5_Hash *hash)
Definition: hash_3d.c:749
int MMG5_bdryPerm(MMG5_pMesh mesh)
Definition: hash_3d.c:2135
API headers for the mmg3d library.
#define MMG3D_LMAX
Definition: libmmg3d.h:58
int MMG3D_update_xtetra(MMG5_pMesh mesh)
Definition: mmg3d2.c:1309
void MMG5_freeXPrisms(MMG5_pMesh mesh)
Definition: zaldy_3d.c:378
void MMG5_freeXTets(MMG5_pMesh mesh)
Definition: zaldy_3d.c:359
static const uint8_t MMG5_idir[4][3]
idir[i]: vertices of face opposite to vertex i
#define MG_REQ
#define MG_GEO
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
#define MMG5_EPSOK
#define MG_EOK(pt)
#define MG_OPNBDY
#define MG_PARBDY
#define MG_MAX(a, b)
#define MMG5_GAP
#define MMG5_ADD_MEM(mesh, size, message, law)
#define MG_EDG(tag)
static const uint8_t MMG5_iprv2[3]
#define MG_EDG_OR_NOM(tag)
#define MMG5_ANGLIM
#define MMG5_TAB_RECALLOC(mesh, ptr, initSize, wantedGap, type, message, law)
#define MMG5_LMAX
static const uint8_t MMG5_inxt2[6]
#define MG_VOK(ppt)
void MMG5_dotprod(int8_t dim, double *a, double *b, double *result)
Definition: tools.c:265
int MMG5_nortri(MMG5_pMesh mesh, MMG5_pTria pt, double *n)
Definition: tools.c:209
#define MG_CRN
double MMG5_orvol(MMG5_pPoint point, MMG5_int *v)
Definition: tools.c:951
#define MG_BDY
#define MMG5_SAFE_MALLOC(ptr, size, type, law)
#define MG_NOSURF
#define MG_NOM
#define MMG5_EPSD2
#define MG_REF
#define MG_PARBDYBDY
#define MMG5_SAFE_FREE(ptr)
#define MMG5_DEL_MEM(mesh, ptr)
MMG5_hgeom * geom
Definition: libmmgtypes.h:575
Identic as MMG5_HGeom but use MMG5_hedge to store edges instead of MMG5_hgeom (memory economy).
Definition: libmmgtypes.h:595
MMG5_hedge * item
Definition: libmmgtypes.h:597
int8_t iso
Definition: libmmgtypes.h:534
int8_t ddebug
Definition: libmmgtypes.h:532
uint8_t nosurf
Definition: libmmgtypes.h:546
double dhd
Definition: libmmgtypes.h:518
int8_t xreg
Definition: libmmgtypes.h:531
int8_t nreg
Definition: libmmgtypes.h:530
MMG mesh structure.
Definition: libmmgtypes.h:605
MMG5_int nc1
Definition: libmmgtypes.h:615
MMG5_Info info
Definition: libmmgtypes.h:651
MMG5_int * adjapr
Definition: libmmgtypes.h:632
MMG5_int ne
Definition: libmmgtypes.h:612
MMG5_pPoint point
Definition: libmmgtypes.h:641
MMG5_int mark
Definition: libmmgtypes.h:618
MMG5_int * adja
Definition: libmmgtypes.h:624
MMG5_int npmax
Definition: libmmgtypes.h:612
MMG5_pxPoint xpoint
Definition: libmmgtypes.h:642
MMG5_HGeom htab
Definition: libmmgtypes.h:650
MMG5_int xp
Definition: libmmgtypes.h:620
MMG5_int base
Definition: libmmgtypes.h:616
MMG5_pTetra tetra
Definition: libmmgtypes.h:643
MMG5_int nt
Definition: libmmgtypes.h:612
MMG5_pTria tria
Definition: libmmgtypes.h:647
MMG5_int np
Definition: libmmgtypes.h:612
MMG5_pxTetra xtetra
Definition: libmmgtypes.h:644
MMG5_int nprism
Definition: libmmgtypes.h:613
MMG5_int na
Definition: libmmgtypes.h:612
MMG5_int * adjt
Definition: libmmgtypes.h:628
MMG5_int xpmax
Definition: libmmgtypes.h:612
Structure to store points of a MMG mesh.
Definition: libmmgtypes.h:270
double n[3]
Definition: libmmgtypes.h:272
int16_t tag
Definition: libmmgtypes.h:284
MMG5_int tmp
Definition: libmmgtypes.h:280
double c[3]
Definition: libmmgtypes.h:271
MMG5_int xp
Definition: libmmgtypes.h:279
MMG5_int s
Definition: libmmgtypes.h:283
MMG5_int flag
Definition: libmmgtypes.h:282
MMG5_int v[4]
Definition: libmmgtypes.h:403
MMG5_int xt
Definition: libmmgtypes.h:407
MMG5_int edg[3]
Definition: libmmgtypes.h:339
int16_t tag[3]
Definition: libmmgtypes.h:342
MMG5_int ref
Definition: libmmgtypes.h:335
MMG5_int flag
Definition: libmmgtypes.h:341
MMG5_int base
Definition: libmmgtypes.h:336
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 ftag[4]
Definition: libmmgtypes.h:423