Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
libmmg3d_tools.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
33#include "libmmg3d.h"
34#include "mmgcommon_private.h"
36#include "mmgversion.h"
38#include "mmgexterns_private.h"
39
40void MMG5_argv_cleanup( char **mmgArgv, int mmgArgc )
41{
42 int i;
43 for ( i = 0; i < mmgArgc; ++i )
44 MMG5_SAFE_FREE(mmgArgv[i]);
45 MMG5_SAFE_FREE(mmgArgv);
46}
47
49
50 if ( mesh->info.ani || (met && met->size == 6) ) {
51 /* Force data consistency: if aniso metric is provided, met->size==6 and
52 * info.ani==0; with -A option, met->size==1 and info.ani==1 */
53 met->size = 6;
54 mesh->info.ani = 1;
55
56 if ( (!met->m) && (!mesh->info.optim) && mesh->info.hsiz<=0. ) {
57 MMG5_caltet = MMG5_caltet_iso;
58 MMG5_caltri = MMG5_caltri_iso;
60 // same as MMG5_lenSurfEdg_iso for iso metric. The edge can be boundary or
61 // intenal but the test relies on the MG_BDY tag that may be missing along
62 // boundary edges (it doesn't matter in iso mode as we always compute the
63 // "straight" edge length). It starts from tetra pointer and edge index.
64 MMG5_lenedg = MMG5_lenedg_iso;
65 // has to be used to compute "straight" edge length from coordinates
67 // Straight edge length (edge is guessed to be a surface edge) from point indices
68 MMG5_lenSurfEdg = MMG5_lenSurfEdg_iso;
69 }
70 else {
71 MMG5_caltet = MMG5_caltet_ani;
72 MMG5_caltri = MMG5_caltri_ani;
74 // lenedg is meant to compute the curve length along boundary edges
75 // and the straight length for internal edges (from
76 // tetra pointer and an edge index) but it relies
77 // on the MG_BDY tag which is not always consistent to detect boundary
78 // edges. Moreover, it seems that the "straight" length is computed in iso
79 // mode - the origin and effect of computing curve lengths along boudary
80 // edges should be investigated...
81 MMG5_lenedg = MMG5_lenedg_ani;
82 // lenedgCoor has to be used to compute "straight" edge length from
83 // coordinates
85 // lenSurfEdg can be called only from a boundary edge: curve length for
86 // aniso metric from point indices
87 MMG5_lenSurfEdg = MMG5_lenSurfEdg_ani;
88 }
89 MMG5_intmet = MMG5_intmet_ani;
90 // warning the lenedg_ani function we may erroneously approximate the length
91 // of a curve boundary edge by the length of the straight edge if the
92 // "MG_BDY" tag is missing along the edge.
93 MMG5_lenedgspl = MMG5_lenedg_ani;
94 MMG5_movintpt = MMG5_movintpt_ani;
95 MMG5_movbdyregpt = MMG5_movbdyregpt_ani;
96 MMG5_movbdyrefpt = MMG5_movbdyrefpt_ani;
97 MMG5_movbdynompt = MMG5_movbdynompt_ani;
98 MMG5_movbdyridpt = MMG5_movbdyridpt_ani;
99 MMG5_interp4bar = MMG5_interp4bar_ani;
100 MMG5_compute_meanMetricAtMarkedPoints = MMG5_compute_meanMetricAtMarkedPoints_ani;
101 MMG3D_defsiz = MMG3D_defsiz_ani;
102 MMG3D_gradsiz = MMG3D_gradsiz_ani;
103 MMG3D_gradsizreq = MMG3D_gradsizreq_ani;
104#ifndef MMG_PATTERN
105 MMG5_cavity = MMG5_cavity_ani;
106 MMG3D_PROctreein = MMG3D_PROctreein_ani;
107#endif
108 }
109 else {
110 if ( mesh->info.optimLES ) {
111 MMG5_caltet = MMG3D_caltetLES_iso;
112 MMG5_movintpt = MMG5_movintpt_iso;
113 }
114 else {
115 MMG5_caltet = MMG5_caltet_iso;
116 MMG5_movintpt = MMG5_movintpt_iso;
117 }
118 MMG5_caltri = MMG5_caltri_iso;
120 // same as MMG5_lenSurfEdg_iso for iso metric. The edge can be boundary or
121 // intenal but the test relies on the MG_BDY tag that may be missing along
122 // boundary edges (it doesn't matter in iso mode as we always compute the
123 // "straight" edge length). It starts from tetra pointer and edge index.
124 MMG5_lenedg = MMG5_lenedg_iso;
125 // lenedgCoor has to be used to compute "straight" edge length from coordinates
127 // Straight edge length (edge is guessed to be a surface edge) from point indices
128 MMG5_lenSurfEdg = MMG5_lenSurfEdg_iso;
129 MMG5_intmet = MMG5_intmet_iso;
130 MMG5_lenedgspl = MMG5_lenedg_iso;
131 MMG5_movbdyregpt = MMG5_movbdyregpt_iso;
132 MMG5_movbdyrefpt = MMG5_movbdyrefpt_iso;
133 MMG5_movbdynompt = MMG5_movbdynompt_iso;
134 MMG5_movbdyridpt = MMG5_movbdyridpt_iso;
135 MMG5_interp4bar = MMG5_interp4bar_iso;
136 MMG5_compute_meanMetricAtMarkedPoints = MMG5_compute_meanMetricAtMarkedPoints_iso;
137 MMG3D_defsiz = MMG3D_defsiz_iso;
138 MMG3D_gradsiz = MMG3D_gradsiz_iso;
139 MMG3D_gradsizreq = MMG3D_gradsizreq_iso;
140
141#ifndef MMG_PATTERN
142 MMG5_cavity = MMG5_cavity_iso;
143 MMG3D_PROctreein = MMG3D_PROctreein_iso;
144#endif
145 }
146}
147
148int MMG3D_Get_adjaTet(MMG5_pMesh mesh, MMG5_int kel, MMG5_int listet[4]) {
149 MMG5_int idx;
150
151 if ( ! mesh->adja ) {
152 if (! MMG3D_hashTetra(mesh, 0))
153 return 0;
154 }
155
156 idx = 4*(kel-1);
157 listet[0] = mesh->adja[idx+1]/4;
158 listet[1] = mesh->adja[idx+2]/4;
159 listet[2] = mesh->adja[idx+3]/4;
160 listet[3] = mesh->adja[idx+4]/4;
161
162 return 1;
163}
164
165int MMG3D_usage(char *prog) {
166
167 /* Common generic options, file options and mode options */
168 MMG5_mmgUsage(prog);
169
170 /* Lagrangian option (only for mmg2d/3d) */
172
173 /* Common parameters (first section) */
175
176 /* Parameters shared by mmg2d and 3d only*/
178
179#ifndef MMG_PATTERN
180 fprintf(stdout,"-octree val specify the max number of points per octree cell \n");
181#endif
182#ifdef USE_SCOTCH
183 fprintf(stdout,"-rn [n] turn on or off the renumbering using SCOTCH [1/0] \n");
184#endif
185 fprintf(stdout,"\n");
186
187 fprintf(stdout,"-nofem do not force Mmg to create a finite element mesh \n");
188 fprintf(stdout,"-nosurf no surface modifications\n");
189
190 fprintf(stdout,"\n");
191
192 /* Common parameters (second section) */
194
195 fprintf(stdout,"-optimLES enable skewness improvement (for LES computations)\n");
196
197 /* Common options for advanced users */
199
200 fprintf(stdout,"\n\n");
201
202 return 1;
203}
204
206
208
209#ifndef MMG_PATTERN
210 fprintf(stdout,"Max number of point per octree cell (-octree) : %d\n",
212#endif
213#ifdef USE_SCOTCH
214 fprintf(stdout,"SCOTCH renumbering : enabled\n");
215#else
216 fprintf(stdout,"SCOTCH renumbering : disabled\n");
217#endif
218 fprintf(stdout,"\n\n");
219
220 return 1;
221}
222
238int MMG3D_storeknownar(int argc,char *argv[],MMG5_pMesh mesh,MMG5_pSol met,
239 MMG5_pSol sol,int *mmgArgc, char *mmgArgv[]) {
240 MMG5_pSol tmp = NULL;
241 double val;
242 int i;
243 char namein[MMG5_FILESTR_LGTH],*endptr;
244 int param;
245
246 i = 1;
247 while ( i < argc ) {
248 if ( *argv[i] == '-' ) {
249 switch(argv[i][1]) {
250 case 'a':
251 if ( !strcmp(argv[i],"-ar") ) {
252 if ( i >= argc -1 ) {
253 fprintf(stderr,"\nMissing argument option %s\n",argv[i]);
254 return 0;
255 }
256 else {
257 val = strtof(argv[i+1],&endptr);
258 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
259 ++i;
261 return 0;
262 }
263 else {
264 /* argument is not a number */
265 fprintf(stderr,"\nMissing argument option %s\n",argv[i]);
266 return 0;
267 }
268 }
269 }
270 else {
271 /* Arg unknown by Mmg: arg starts with -a but is not known */
272 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
273 }
274 break;
275 case 'A': /* anisotropy */
276 if ( !strcmp(argv[i],"-A") ) {
278 return 0;
279 }
280 else {
281 /* Arg unknown by Mmg: arg starts with -A but is not known */
282 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
283 }
284 break;
285 case 'd':
286 if ( !strcmp(argv[i],"-default") ) {
287 mesh->mark=1;
288 }
289 else if ( !strcmp(argv[i],"-d") ) {
290 /* debug */
292 return 0;
293 }
294 }
295 else {
296 /* Arg unknown by Mmg: arg starts with -d but is not known */
297 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
298 }
299 break;
300 case 'f':
301 if ( !strcmp(argv[i],"-f") ) {
302 if ( ++i < argc && isascii(argv[i][0]) && argv[i][0]!='-' ) {
303 if ( !MMG3D_Set_inputParamName(mesh,argv[i]) )
304 return 0;
305 }
306 else {
307 fprintf(stderr,"\nMissing filename for %s\n",argv[i-1]);
308 return 0;
309 }
310 }
311 else {
312 /* Arg unknown by Mmg: arg starts with -f but is not known */
313 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
314 }
315 break;
316 case 'h':
317 param = MMG5_UNSET;
318 if ( i >= argc -1 ) {
319 fprintf(stderr,"\nMissing argument option %s\n",argv[i]);
320 return 0;
321 }
322 else {
323 if ( !strcmp(argv[i],"-hmin") ) {
324 param = MMG3D_DPARAM_hmin;
325 }
326 else if ( !strcmp(argv[i],"-hmax") ) {
327 param = MMG3D_DPARAM_hmax;
328 }
329 else if ( !strcmp(argv[i],"-hsiz") ) {
330 param = MMG3D_DPARAM_hsiz;
331 }
332 else if ( !strcmp(argv[i],"-hausd") ) {
333 param = MMG3D_DPARAM_hausd;
334 }
335 else if ( !strcmp(argv[i],"-hgradreq") ) {
336 param = MMG3D_DPARAM_hgradreq;
337 }
338 else if ( !strcmp(argv[i],"-hgrad") ) {
339 param = MMG3D_DPARAM_hgrad;
340 }
341 else {
342 /* Arg unknown by Mmg: arg starts with -h but is not known */
343 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
344 }
345
346 if ( param != MMG5_UNSET ) {
347
348 /* Arg can be parsed */
349 val = strtof(argv[i+1],&endptr);
350 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
351 ++i;
352 if ( !MMG3D_Set_dparameter(mesh,met,param,val) ){
353 return 0;
354 }
355 } else {
356 /* argument is not a number */
357 fprintf(stderr,"\nMissing argument option %s\n",argv[i]);
358 return 0;
359 }
360 }
361 }
362
363 break;
364 case 'i':
365 if ( !strcmp(argv[i],"-in") ) {
366 if ( ++i < argc && isascii(argv[i][0]) && argv[i][0]!='-') {
367 if ( !MMG3D_Set_inputMeshName(mesh, argv[i]) )
368 return 0;
369
370 }else{
371 fprintf(stderr,"\nMissing filname for %s\n",argv[i-1]);
372 return 0;
373 }
374 }
375 else if ( !strcmp(argv[i],"-isoref") && ++i <= argc ) {
377 atoi(argv[i])) )
378 return 0;
379 }
380 else {
381 /* Arg unknown by Mmg: arg starts with -i but is not known */
382 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
383 }
384 break;
385 case 'l':
386 if ( !strcmp(argv[i],"-lag") ) {
387 if ( ++i < argc && isdigit(argv[i][0]) ) {
388 if ( !MMG3D_Set_iparameter(mesh,met,MMG3D_IPARAM_lag,atoi(argv[i])) )
389 return 0;
390 }
391 else {
392 fprintf(stderr,"\nMissing or unexpected argument option %s\n",argv[i-1]);
393 MMG3D_usage(argv[0]);
394 return 0;
395 }
396 }
397 else if ( !strcmp(argv[i],"-ls") ) {
399 return 0;
400
401 if ( i < argc -1 ) {
402 val = strtof(argv[i+1],&endptr);
403 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
404 ++i;
406 return 0;
407 }
408 }
409 }
410 else if ( !strcmp(argv[i],"-lssurf") ) {
412 return 0;
413
414 if ( i < argc -1 ) {
415 val = strtof(argv[i+1],&endptr);
416 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
417 ++i;
419 return 0;
420 }
421 }
422 }
423 else {
424 /* Arg unknown by Mmg: arg starts with -l but is not known */
425 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
426 }
427 break;
428 case 'm':
429 if ( !strcmp(argv[i],"-met") ) {
430 if ( !met ) {
431 fprintf(stderr,"\nNo metric structure allocated for %s option\n",
432 argv[i-1]);
433 return 0;
434 }
435 if ( ++i < argc && isascii(argv[i][0]) && argv[i][0]!='-' ) {
436 if ( !MMG3D_Set_inputSolName(mesh,met,argv[i]) )
437 return 0;
438 }
439 else {
440 fprintf(stderr,"\nMissing filname for %s\n",argv[i-1]);
441 return 0;
442 }
443 }
444 else if ( !strcmp(argv[i],"-m") ) {
445 /* memory */
446 if ( ++i < argc && isdigit(argv[i][0]) ) {
447 if ( !MMG3D_Set_iparameter(mesh,met,MMG3D_IPARAM_mem,atoi(argv[i])) )
448 return 0;
449 }
450 else {
451 fprintf(stderr,"\nMissing argument option %s\n",argv[i-1]);
452 return 0;
453 }
454 }
455 else {
456 /* Arg unknown by Mmg: arg starts with -m but is not known */
457 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
458 }
459 break;
460 case 'n':
461 if ( !strcmp(argv[i],"-nofem") ) {
463 return 0;
464 }
465 else if ( !strcmp(argv[i],"-nreg") ) {
467 return 0;
468 }
469 else if ( !strcmp(argv[i],"-nr") ) {
471 return 0;
472 }
473 else if ( !strcmp(argv[i],"-nsd") ) {
474 if ( ++i < argc && isdigit(argv[i][0]) ) {
475 if ( !MMG3D_Set_iparameter(mesh,met,MMG3D_IPARAM_numsubdomain,atoi(argv[i])) )
476 return 0;
477 }
478 else {
479 fprintf(stderr,"\nMissing argument option %s\n",argv[i-1]);
480 return 0;
481 }
482 }
483 else if ( !strcmp(argv[i],"-noswap") ) {
485 return 0;
486 }
487 else if( !strcmp(argv[i],"-noinsert") ) {
489 return 0;
490 }
491 else if( !strcmp(argv[i],"-nomove") ) {
493 return 0;
494 }
495 else if( !strcmp(argv[i],"-nosurf") ) {
497 return 0;
498 }
499 else if( !strcmp(argv[i],"-nosizreq") ) {
501 return 0;
502 }
503 }
504 else {
505 /* Arg unknown by Mmg: arg starts with -n but is not known */
506 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
507 }
508 break;
509 case 'o':
510 if ( (!strcmp(argv[i],"-out")) || (!strcmp(argv[i],"-o")) ) {
511 if ( ++i < argc && isascii(argv[i][0]) && argv[i][0]!='-') {
512 if ( !MMG3D_Set_outputMeshName(mesh,argv[i]) )
513 return 0;
514 }else{
515 fprintf(stderr,"\nMissing filname for %s\n",argv[i-1]);
516 return 0;
517 }
518 }
519 else if ( !strcmp(argv[i],"-opnbdy") ) {
521 return 0;
522 }
523#ifndef MMG_PATTERN
524 else if ( !strcmp(argv[i],"-octree") && ++i < argc ) {
526 atoi(argv[i])) )
527 return 0;
528 }
529#endif
530 else if( !strcmp(argv[i],"-optimLES") ) {
532 return 0;
533 }
534 else if( !strcmp(argv[i],"-optim") ) {
536 return 0;
537 }
538 else {
539 /* Arg unknown by Mmg: arg starts with -o but is not known */
540 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
541 }
542 break;
543 case 'r':
544 if ( !strcmp(argv[i],"-rmc") ) {
546 return 0;
547 if ( i < argc -1 ) {
548 val = strtof(argv[i+1],&endptr);
549 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
550 ++i;
552 return 0;
553 }
554 }
555 }
556#ifdef USE_SCOTCH
557 else if ( !strcmp(argv[i],"-rn") ) {
558 if ( ++i < argc ) {
559 if ( isdigit(argv[i][0]) ) {
560 if ( !MMG3D_Set_iparameter(mesh,met,MMG3D_IPARAM_renum,atoi(argv[i])) )
561 return 0;
562 }
563 else {
564 fprintf(stderr,"\nMissing argument option %s\n",argv[i-1]);
565 return 0;
566 }
567 }
568 else {
569 fprintf(stderr,"\nMissing argument option %s\n",argv[i-1]);
570 return 0;
571 }
572 }
573#endif
574 else {
575 /* Arg unknown by Mmg: arg starts with -r but is not known */
576 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
577 }
578
579 break;
580 case 's':
581 if ( !strcmp(argv[i],"-sol") ) {
582 /* For retrocompatibility, store the metric if no sol structure available */
583 tmp = sol ? sol : met;
584 assert(tmp);
585 if ( ++i < argc && isascii(argv[i][0]) && argv[i][0]!='-' ) {
586 if ( !MMG3D_Set_inputSolName(mesh,tmp,argv[i]) )
587 return 0;
588 }
589 else {
590 fprintf(stderr,"\nMissing filname for %s\n",argv[i-1]);
591 return 0;
592 }
593 }
594 else {
595 /* Arg unknown by Mmg: arg starts with -s but is not known */
596 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
597 }
598 break;
599 case 'v':
600 if ( !strcmp(argv[i],"-v") ) {
601 if ( ++i < argc ) {
602 if ( isdigit(argv[i][0]) ||
603 (argv[i][0]=='-' && isdigit(argv[i][1])) ) {
604 if ( !MMG3D_Set_iparameter(mesh,met,MMG3D_IPARAM_verbose,atoi(argv[i])) )
605 return 0;
606 }
607 else {
608 i--;
609 fprintf(stderr,"\nMissing argument option %s\n",argv[i]);
610 }
611 }
612 else {
613 fprintf(stderr,"\nMissing argument option %s\n",argv[i-1]);
614 return 0;
615 }
616 }
617 else {
618 /* Arg unknown by Mmg: arg starts with -v but is not known */
619 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
620 }
621 break;
622 case 'x':
623 if ( !strcmp(argv[i],"-xreg") ) {
625 return 0;
626 if ( i < argc -1 ) {
627 val = strtof(argv[i+1],&endptr);
628 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
629 ++i;
631 return 0;
632 }
633 }
634 }
635 else {
636 /* Arg unknown by Mmg: arg starts with -x but is not known */
637 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
638 }
639 break;
640 default:
641 /* Arg unknown by Mmg: arg starts with -<.>, <.> being a letter with
642 * which no known argument of Mmg begins */
643 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
644 }
645 }
646 else {
647 /* Arg unknown by Mmg: arg doesn't start by "-" and has not been parsed as
648 * the value expected by another argument */
649 MMG_ARGV_APPEND(argv, mmgArgv, i, *mmgArgc,return 0);
650 }
651 i++;
652 }
653
654 return 1;
655}
656
657
658// In ls mode : metric must be provided using -met option (-sol or default is the ls).
659// In adp mode : -sol or -met or default allow to store the metric.
660int MMG3D_parsar(int argc,char *argv[],MMG5_pMesh mesh,MMG5_pSol met,MMG5_pSol sol) {
661 MMG5_pSol tmp = NULL;
662 int i;
663 char namein[MMG5_FILESTR_LGTH];
664
665 /* First step: search if user want to see the default parameters values or is
666 * asking for help */
667 for ( i=1; i< argc; ++i ) {
668 if ( !strcmp(argv[i],"-val") ) {
669 if ( !MMG3D_defaultValues(mesh) ) return 0;
670 return 0;
671 }
672 else if ( ( !strcmp( argv[ i ],"-?" ) ) || ( !strcmp( argv[ i ],"-h" ) ) ) {
673 MMG3D_usage(argv[0]);
674 return 0;
675 }
676 }
677
678 /* Second step: read all other arguments known by Mmg and append unknown ones
679 * in mmgArgv. */
680 int mmgArgc = 0;
681 char** mmgArgv = NULL;
682 MMG5_SAFE_MALLOC( mmgArgv, argc, char*,return 0);
683 MMG_ARGV_APPEND ( argv, mmgArgv, 0, mmgArgc,return 0);
684
685 int ier = MMG3D_storeknownar( argc,argv,mesh,met,sol,&mmgArgc,mmgArgv);
686
687 /* Third step: treat unknown args */
688 if ( ier ) {
689 i = 1;
690 while ( i < mmgArgc ) {
691 if ( *mmgArgv[i] != '-' ) {
692 /* Arg doesn't start by '-', try to parse it as filname */
693 if ( mesh->namein == NULL ) {
694 if ( !MMG3D_Set_inputMeshName(mesh,mmgArgv[i]) )
695 return 0;
696 }
697 else if ( mesh->nameout == NULL ) {
698 if ( !MMG3D_Set_outputMeshName(mesh,mmgArgv[i]) )
699 return 0;
700 }
701 else {
702 fprintf(stdout,"Argument %s ignored\n",mmgArgv[i]);
703 ier = 0;
704 break;
705 }
706 }
707 else {
708 /* Arg start by '-' and has not been parsed: unexpected */
709 fprintf(stdout,"Argument %s ignored\n",mmgArgv[i]);
710 ier = 0;
711 break;
712 }
713 i++;
714 }
715 }
716 MMG5_argv_cleanup(mmgArgv,mmgArgc);
717
718 /* Fourth step: handle errors */
719 if ( !ier ) {
720 MMG3D_usage(argv[0]);
721 return 0;
722 }
723
724 /* Last step: check file names */
725 if ( mesh->info.imprim == -99 ) {
726 fprintf(stdout,"\n -- PRINT (0 10(advised) -10) ?\n");
727 fflush(stdin);
728 MMG_FSCANF(stdin,"%d",&i);
730 return 0;
731 }
732
733 if ( mesh->namein == NULL ) {
734 fprintf(stdout,"\n -- INPUT MESH NAME ?\n");
735 fflush(stdin);
736 MMG_FSCANF(stdin,"%127s",namein);
737 if ( !MMG3D_Set_inputMeshName(mesh,namein) )
738 return 0;
739 }
740
741 if ( mesh->nameout == NULL ) {
743 return 0;
744 }
745
746 /* adp mode: if the metric name has been stored in sol, move it in met */
747 if ( met->namein==NULL && sol && sol->namein && !(mesh->info.iso || mesh->info.isosurf || mesh->info.lag>=0) ) {
749 return 0;
751 }
752
753 /* default : store solution (resp. displacement) name in iso
754 * (resp. lagrangian) mode, metric name otherwise */
755 tmp = ( mesh->info.iso || mesh->info.isosurf || mesh->info.lag >=0 ) ? sol : met;
756 assert ( tmp );
757 if ( tmp->namein == NULL ) {
758 if ( !MMG3D_Set_inputSolName(mesh,tmp,"") ) { return 0; }
759 }
760 if ( met->nameout == NULL ) {
761 if ( !MMG3D_Set_outputSolName(mesh,met,"") )
762 return 0;
763 }
764
765 return 1;
766}
767
769 float fp1,fp2,hausd;
770 int i,j,ret,npar,nbr,split;
771 MMG5_int ref,rin,rex,br;
772 char *ptr,buf[256],data[256];
773 FILE *in;
774 fpos_t position;
775
776 /* check for parameter file */
777 if (mesh->info.fparam) {
778 strcpy(data,mesh->info.fparam);
779 }
780 else {
781 strcpy(data,mesh->namein);
782 }
783
784 ptr = MMG5_Get_filenameExt(data);
785
786 if ( ptr ) *ptr = '\0';
787 strcat(data,".mmg3d");
788
789 in = fopen(data,"rb");
790 if ( !in ) {
791 if ( !mesh->info.fparam ) {
792 strcat(data,".mmg3d5");
793 in = fopen(data,"rb");
794 if ( !in ) {
795 sprintf(data,"%s","DEFAULT.mmg3d");
796 in = fopen(data,"rb");
797 if ( !in ) {
798 sprintf(data,"%s","DEFAULT.mmg3d5");
799 in = fopen(data,"rb");
800 if ( !in ) {
801 return 1;
802 }
803 }
804 }
805 }
806 else if (mesh->info.fparam ) {
807 fprintf(stderr," ** In %s: %s file NOT FOUND. \n",__func__,data);
808 fprintf(stdout," ## ERROR: UNABLE TO LOAD PARAMETER FILE.\n");
809 return 0;
810 }
811 }
812 if ( mesh->info.imprim >= 0 )
813 fprintf(stdout,"\n %%%% %s OPENED\n",data);
814
815 /* read parameters */
816 while ( !feof(in) ) {
817 /* scan line */
818 ret = fscanf(in,"%255s",data);
819 if ( !ret || feof(in) ) break;
820 for (i=0; i<strlen(data); i++) data[i] = tolower(data[i]);
821
822 /* Read user-defined references for the LS mode */
823 if ( !strcmp(data,"lsreferences") ) {
824 ret = fscanf(in,"%d",&npar);
825 if ( !ret ) {
826 fprintf(stderr," %%%% Wrong format for lsreferences: %d\n",npar);
827 return 0;
828 }
829
831 return 0;
832 }
833 for (i=0; i<mesh->info.nmat; i++) {
834 MMG_FSCANF(in,"%" MMG5_PRId "",&ref);
835 fgetpos(in,&position);
836 MMG_FSCANF(in,"%255s",data);
837 split = MMG5_MMAT_NoSplit;
838 rin = rex = ref;
839 if ( strcmp(data,"nosplit") ) {
840 fsetpos(in,&position);
841 split = MMG5_MMAT_Split;
842 MMG_FSCANF(in,"%" MMG5_PRId "",&rin);
843 MMG_FSCANF(in,"%" MMG5_PRId "",&rex);
844 }
845 if ( !MMG3D_Set_multiMat(mesh,met,ref,split,rin,rex) ) {
846 return 0;
847 }
848 }
849 }
850 /* Read user-defined local parameters and store them in the structure info->par */
851 else if ( !strcmp(data,"parameters") ) {
852 ret = fscanf(in,"%d",&npar);
853
854 if ( !ret ) {
855 fprintf(stderr," %%%% Wrong format for parameters: %d\n",npar);
856 return 0;
857 }
858
860 return 0;
861 }
862
863 for (i=0; i<mesh->info.npar; i++) {
864 ret = fscanf(in,"%" MMG5_PRId " %255s ",&ref,buf);
865 if ( ret )
866 ret = fscanf(in,"%f %f %f",&fp1,&fp2,&hausd);
867
868 if ( !ret ) {
869 fprintf(stderr," %%%% Wrong format: %s\n",buf);
870 return 0;
871 }
872
873 for (j=0; j<strlen(buf); j++) buf[j] = tolower(buf[j]);
874
875 if ( (!strcmp(buf,"triangles") || !strcmp(buf,"triangle")) ) {
876 if ( !MMG3D_Set_localParameter(mesh,met,MMG5_Triangle,ref,fp1,fp2,hausd) ) {
877 return 0;
878 }
879 }
880 else if ( !strcmp(buf,"tetrahedra") || !strcmp(buf,"tetrahedron") ) {
881 if ( !MMG3D_Set_localParameter(mesh,met,MMG5_Tetrahedron,ref,fp1,fp2,hausd) ) {
882 return 0;
883 }
884 }
885 else {
886 fprintf(stderr," %%%% Wrong format: %s\n",buf);
887 return 0;
888 }
889 }
890 }
891 else if ( !strcmp(data,"lsbasereferences") ) {
892 MMG_FSCANF(in,"%d",&nbr);
894 return 0;
895
896 for (i=0; i<mesh->info.nbr; i++) {
897 MMG_FSCANF(in,"%" MMG5_PRId "",&br);
898 if ( !MMG3D_Set_lsBaseReference(mesh,met,br) ) {
899 return 0;
900 }
901 }
902 }
903 }
904 fclose(in);
905 return 1;
906}
907
909
910 free(mesh->info.par);
911 mesh->info.npar = 0;
912
913 return 1;
914}
915
917 MMG5_pTetra pt,pt1;
918 MMG5_pPrism pp;
919 MMG5_pTria ptt;
920 MMG5_Hash hash;
921 int i;
922 MMG5_int ref,*adja,j,k,iel;
923
924 *nb_tria = 0;
925 memset ( &hash, 0x0, sizeof(MMG5_Hash));
926
927 if ( !mesh->tetra ) {
928 /* No triangle at all */
929 return 1;
930 }
931
934 if ( !mesh->adja ) {
935 /* create tetra adjacency */
936 if ( !MMG3D_hashTetra( mesh,0 ) ) {
937 fprintf(stderr,"\n ## Error: %s: unable to create "
938 "adjacency table.\n",__func__);
939 return 0;
940 }
941 }
942
943 if ( !mesh->adjapr ) {
944 /* create prism adjacency */
945 if ( !MMG3D_hashPrism(mesh) ) {
946 fprintf(stderr,"\n ## Error: %s: Prism hashing problem.\n",__func__);
947 return 0;
948 }
949 }
950
951 /* If mesh->xtetra is filled, we assume that the surface analysis is
952 * complete */
953 if ( !mesh->xtetra ) {
954 /* compatibility triangle orientation w/r tetras */
955 if ( !MMG5_bdryPerm(mesh) ) {
956 fprintf(stderr,"\n ## Error: %s: Boundary orientation problem.\n",__func__);
957 return 0;
958 }
959 /* identify surface mesh */
960 if ( !MMG5_chkBdryTria(mesh) ) {
961 fprintf(stderr,"\n ## Error: %s: Boundary problem.\n",__func__);
962 return 0;
963 }
966
967 /* create surface adjacency */
968 if ( !MMG3D_hashTria(mesh,&hash) ) {
969 MMG5_DEL_MEM(mesh,hash.item);
970 fprintf(stderr,"\n ## Error: %s: Hashing problem.\n",__func__);
971 return 0;
972 }
973
974 /* identify connexity and flip orientation of faces if needed */
975 if ( !MMG5_setadj(mesh) ) {
976 fprintf(stderr,"\n ## Error: %s: Topology problem.\n",__func__);
977 MMG5_DEL_MEM(mesh,hash.item);
978 return 0;
979 }
980
981 /* set bdry entities to tetra and fill the orientation field */
982 if ( !MMG5_bdrySet(mesh) ) {
983 MMG5_DEL_MEM(mesh,hash.item);
984 fprintf(stderr,"\n ## Error: %s: Boundary problem.\n",__func__);
985 return 0;
986 }
987 MMG5_DEL_MEM(mesh,hash.item);
988 }
989
991 for ( k=1; k<=mesh->ne; k++ ) {
992 pt = &mesh->tetra[k];
993 if ( !MG_EOK(pt) ) continue;
994
995 adja = &mesh->adja[4*(k-1)+1];
996
997 for ( i=0; i<4; i++ ) {
998 iel = adja[i] / 4;
999 assert ( iel != k );
1000
1001 pt1 = &mesh->tetra[iel];
1002
1003 if ( (!iel) || (pt->ref != pt1->ref) ||
1004 (mesh->info.opnbdy && pt->xt &&
1005 (mesh->xtetra[pt->xt].ftag[i] & MG_BDY) ) ) {
1006 /* Do not treat boundary faces */
1007 continue;
1008 }
1009 if ( k < iel ) {
1010 /* Treat face from the tetra with lowest index */
1011 ++(*nb_tria);
1012 }
1013 }
1014 }
1015 for ( k=1; k<=mesh->nprism; k++ ) {
1016 pp = &mesh->prism[k];
1017 if ( !MG_EOK(pp) ) continue;
1018
1019 adja = &mesh->adjapr[5*(k-1)+1];
1020
1021 for ( i=0; i<2; i++ ) {
1022 iel = adja[i] / 5;
1023
1024 if ( iel<0 ) {
1025 ref = mesh->tetra[MMG5_abs(iel)].ref;
1026 } else {
1027 ref = mesh->prism[iel].ref;
1028 }
1029
1030 if ( (!iel) || (pp->ref != ref) ||
1031 (mesh->info.opnbdy && pp->xpr &&
1032 (mesh->xprism[pp->xpr].ftag[i] & MG_BDY) ) ) {
1033 /* Do not treat boundary faces */
1034 continue;
1035 }
1036 if ( k < iel ) {
1037 /* Treat face from the element with lowest index */
1038 ++(*nb_tria);
1039 }
1040 }
1041 }
1042
1043 if ( !(*nb_tria) ) {
1044 return 1;
1045 }
1046
1048 if ( mesh->nt ) {
1049 MMG5_ADD_MEM(mesh,(*nb_tria)*sizeof(MMG5_Tria),"non boundary triangles",
1050 printf(" Exit program.\n");
1051 MMG5_DEL_MEM(mesh,hash.item);
1052 return 0);
1053 MMG5_SAFE_RECALLOC(mesh->tria,(mesh->nt+1),(mesh->nt+(*nb_tria)+1),
1054 MMG5_Tria,"non bdy tria arrray",return 0);
1055 }
1056 else {
1057 MMG5_ADD_MEM(mesh,((*nb_tria)+1)*sizeof(MMG5_Tria),"non boundary triangles",
1058 printf(" Exit program.\n");
1059 MMG5_DEL_MEM(mesh,hash.item);
1060 return 0);
1061 MMG5_SAFE_RECALLOC(mesh->tria,0,((*nb_tria)+1),
1062 MMG5_Tria,"non bdy tria arrray",return 0);
1063 }
1064
1065 j = mesh->nt+1;
1066 for ( k=1; k<=mesh->ne; k++ ) {
1067 pt = &mesh->tetra[k];
1068 if ( !MG_EOK(pt) ) continue;
1069
1070 adja = &mesh->adja[4*(k-1)+1];
1071
1072 for ( i=0; i<4; i++ ) {
1073 iel = adja[i] / 4;
1074 assert ( iel != k );
1075
1076 pt1 = &mesh->tetra[iel];
1077
1078 if ( (!iel) || (pt->ref != pt1->ref) ||
1079 (mesh->info.opnbdy && pt->xt &&
1080 (mesh->xtetra[pt->xt].ftag[i] & MG_BDY)) ) {
1081 /* Do not treat boundary faces */
1082 continue;
1083 }
1084 if ( k < iel ) {
1085 /* Treat edge from the triangle with lowest index */
1086 ptt = &mesh->tria[j++];
1087 assert ( ptt );
1088 ptt->v[0] = pt->v[MMG5_idir[i][0]];
1089 ptt->v[1] = pt->v[MMG5_idir[i][1]];
1090 ptt->v[2] = pt->v[MMG5_idir[i][2]];
1091 ptt->ref = mesh->xtetra[pt->xt].ref[i];
1092 }
1093 }
1094 }
1095
1096 for ( k=1; k<=mesh->nprism; k++ ) {
1097 pp = &mesh->prism[k];
1098 if ( !MG_EOK(pp) ) continue;
1099
1100 adja = &mesh->adjapr[5*(k-1)+1];
1101
1102 for ( i=0; i<2; i++ ) {
1103 iel = adja[i] / 5;
1104
1105 if ( iel<0 ) {
1106 ref = mesh->tetra[MMG5_abs(iel)].ref;
1107 } else {
1108 ref = mesh->prism[iel].ref;
1109 }
1110 if ( (!iel) || (pp->ref != ref) ||
1111 (mesh->info.opnbdy && pp->xpr &&
1112 (mesh->xprism[pp->xpr].ftag[i] & MG_BDY)) ) {
1113 /* Do not treat boundary faces */
1114 continue;
1115 }
1116 if ( k < iel ) {
1117 /* Treat edge from the triangle with lowest index */
1118 ptt = &mesh->tria[j++];
1119 assert ( ptt );
1120 ptt->v[0] = pp->v[MMG5_idir_pr[i][0]];
1121 ptt->v[1] = pp->v[MMG5_idir_pr[i][1]];
1122 ptt->v[2] = pp->v[MMG5_idir_pr[i][2]];
1123 ptt->ref = mesh->xprism[pp->xpr].ref[i];
1124 }
1125 }
1126 }
1127
1128 return 1;
1129}
1130
1131int MMG3D_Get_nonBdyTriangle(MMG5_pMesh mesh,MMG5_int* v0,MMG5_int* v1,MMG5_int* v2,
1132 MMG5_int* ref,MMG5_int idx) {
1133 MMG5_pTria ptt;
1134 size_t nt_tot=0;
1135 char *ptr_c = (char*)mesh->tria;
1136
1137 if ( !mesh->tria ) {
1138 fprintf(stderr,"\n ## Error: %s: triangle array is not allocated.\n"
1139 " Please, call the MMG3D_Get_numberOfNonBdyTriangles function"
1140 " before the %s one.\n",
1141 __func__,__func__);
1142 return 0;
1143 }
1144
1145 ptr_c = ptr_c-sizeof(size_t);
1146 nt_tot = *((size_t*)ptr_c);
1147
1148 if ( mesh->nt==(MMG5_int)nt_tot ) {
1149 fprintf(stderr,"\n ## Error: %s: no internal triangle.\n"
1150 " Please, call the MMG3D_Get_numberOfNonBdyTriangles function"
1151 " before the %s one and check that the number of internal"
1152 " triangles is non null.\n",
1153 __func__,__func__);
1154 return 0;
1155 }
1156
1157 if ( mesh->nt+idx > (MMG5_int)nt_tot ) {
1158 fprintf(stderr,"\n ## Error: %s: Can't get the internal triangle of index %" MMG5_PRId "."
1159 " Index must be between 1 and %zu.\n",
1160 __func__,idx,nt_tot-(size_t)mesh->nt);
1161 return 0;
1162 }
1163
1164 ptt = &mesh->tria[mesh->nt+idx];
1165
1166 *v0 = ptt->v[0];
1167 *v1 = ptt->v[1];
1168 *v2 = ptt->v[2];
1169
1170 if ( ref != NULL ) {
1171 *ref = ptt->ref;
1172 }
1173
1174 return 1;
1175}
1176
1178
1179 memcpy(&mesh->info,info,sizeof(MMG5_Info));
1181 if( mesh->info.mem > 0) {
1182 if((mesh->npmax < mesh->np || mesh->ntmax < mesh->nt || mesh->nemax < mesh->ne)) {
1183 return 0;
1184 } else if(mesh->info.mem < 39)
1185 return 0;
1186 }
1187 return 1;
1188}
1189
1191
1192 memcpy(info,&mesh->info,sizeof(MMG5_Info));
1193 return;
1194}
1195
1196int MMG3D_mmg3dcheck(MMG5_pMesh mesh,MMG5_pSol met,MMG5_pSol sol,double critmin, double lmin,
1197 double lmax, MMG5_int *eltab,int8_t metRidTyp) {
1198
1199 mytime ctim[TIMEMAX];
1200 int ier;
1201 char stim[32];
1202
1204
1208
1209 signal(SIGABRT,MMG5_excfun);
1210 signal(SIGFPE,MMG5_excfun);
1211 signal(SIGILL,MMG5_excfun);
1212 signal(SIGSEGV,MMG5_excfun);
1213 signal(SIGTERM,MMG5_excfun);
1214 signal(SIGINT,MMG5_excfun);
1215
1216 tminit(ctim,TIMEMAX);
1217 chrono(ON,&(ctim[0]));
1218
1219 /* Check options */
1220 if ( mesh->info.lag > -1 ) {
1221 fprintf(stderr,"\n ## Error: lagrangian mode unavailable (MMG3D_IPARAM_lag):\n"
1222 " You must call the MMG3D_mmg3dmov function to move a rigidbody.\n");
1224 }
1225 else if ( mesh->info.iso ) {
1226 fprintf(stderr,"\n ## Error: level-set discretisation unavailable"
1227 " (MMG3D_IPARAM_iso):\n"
1228 " You must call the MMG3D_mmg3dmov function to use this option.\n");
1230 }
1231
1232#ifdef USE_SCOTCH
1234#endif
1235
1236 if ( mesh->info.imprim > 0 ) fprintf(stdout,"\n -- MMG3DLIB: INPUT DATA\n");
1237 /* load data */
1238 chrono(ON,&(ctim[1]));
1240
1241 if ( met ) {
1242 if ( met->np && (met->np != mesh->np) ) {
1243 fprintf(stdout," ## WARNING: WRONG SOLUTION NUMBER. IGNORED\n");
1244 MMG5_DEL_MEM(mesh,met->m);
1245 met->np = 0;
1246 }
1247 else if ( met->size!=1 && met->size!=6 ) {
1248 fprintf(stderr,"\n ## ERROR: WRONG DATA TYPE.\n");
1250 }
1251 }
1252 if ( sol ) {
1253 if ( sol->np && (sol->np != mesh->np) ) {
1254 fprintf(stdout," ## WARNING: WRONG SOLUTION NUMBER. IGNORED\n");
1255 MMG5_DEL_MEM(mesh,sol->m);
1256 sol->np = 0;
1257 }
1258 else if ( sol->size!=1 && sol->size!=6 ) {
1259 fprintf(stderr,"\n ## ERROR: WRONG DATA TYPE.\n");
1261 }
1262 }
1263
1264 chrono(OFF,&(ctim[1]));
1265 printim(ctim[1].gdif,stim);
1266 if ( mesh->info.imprim > 0 )
1267 fprintf(stdout," -- INPUT DATA COMPLETED. %s\n",stim);
1268
1269 /* analysis */
1270 chrono(ON,&(ctim[2]));
1271 if ( met ) {
1272 MMG3D_setfunc(mesh,met);
1273 }
1274
1275 if ( mesh->info.imprim > 0 ) {
1276 fprintf(stdout,"\n %s\n MODULE MMG3D: IMB-LJLL : %s (%s)\n %s\n",
1277 MG_STR,MMG_VERSION_RELEASE,MMG_RELEASE_DATE,MG_STR);
1278 fprintf(stdout,"\n -- PHASE 1 : ANALYSIS\n");
1279 }
1280
1281 /* scaling mesh */
1283
1284
1285 if ( met ) {
1286 MMG3D_searchqua(mesh,met,critmin,eltab,metRidTyp);
1287 ier = MMG3D_searchlen(mesh,met,lmin,lmax,eltab,metRidTyp);
1288 if ( !ier )
1290 }
1291
1293}
1294
1295void MMG3D_searchqua(MMG5_pMesh mesh,MMG5_pSol met,double critmin, MMG5_int *eltab,
1296 int8_t metRidTyp) {
1297 MMG5_pTetra pt;
1298 double rap;
1299 MMG5_int k;
1300
1301 assert ( met );
1302
1303 for (k=1; k<=mesh->ne; k++) {
1304 pt = &mesh->tetra[k];
1305
1306 if( !MG_EOK(pt) )
1307 continue;
1308
1309 if ( (!metRidTyp) && met->m && met->size>1 ) {
1310 rap = MMG3D_ALPHAD * MMG5_caltet33_ani(mesh,met,pt);
1311 }
1312 else {
1313 rap = MMG3D_ALPHAD * MMG5_caltet(mesh,met,pt);
1314 }
1315
1316 if ( rap == 0.0 || rap < critmin ) {
1317 eltab[k] = 1;
1318 }
1319 }
1320 return;
1321}
1322
1323int MMG3D_Get_tetFromTria(MMG5_pMesh mesh, MMG5_int ktri, MMG5_int *ktet, int *iface) {
1324 MMG5_int val;
1325
1326 val = mesh->tria[ktri].cc;
1327
1328 if ( !val ) {
1329 fprintf(stderr," ## Error: %s: the main fonction of the Mmg library must be"
1330 " called before this function.\n",__func__);
1331 return 0;
1332 }
1333
1334 *ktet = val/4;
1335
1336 *iface = val%4;
1337
1338 return 1;
1339}
1340
1341
1342int MMG3D_Get_tetsFromTria(MMG5_pMesh mesh, MMG5_int ktri, MMG5_int ktet[2], int iface[2])
1343{
1344 int ier;
1345 MMG5_int itet;
1346#ifndef NDEBUG
1347 MMG5_int ia0,ib0,ic0,ia1,ib1,ic1;
1348#endif
1349
1350 ktet[0] = ktet[1] = 0;
1351 iface[0] = iface[1] = 0;
1352
1353 ier = MMG3D_Get_tetFromTria(mesh,ktri,ktet,iface);
1354
1355 if ( !ier ) return 0;
1356
1357 if ( !mesh->adja ) {
1358 if (!MMG3D_hashTetra(mesh, 0) )
1359 return 0;
1360 }
1361
1362 itet = mesh->adja[4*(*ktet-1) + *iface + 1 ];
1363
1364 if ( itet ) {
1365 ktet[1] = itet/4;
1366 iface[1] = itet%4;
1367
1368#ifndef NDEBUG
1369 ia0 = mesh->tetra[ktet[0]].v[MMG5_idir[iface[0]][0]];
1370 ib0 = mesh->tetra[ktet[0]].v[MMG5_idir[iface[0]][1]];
1371 ic0 = mesh->tetra[ktet[0]].v[MMG5_idir[iface[0]][2]];
1372
1373 ia1 = mesh->tetra[ktet[1]].v[MMG5_idir[iface[1]][0]];
1374 ib1 = mesh->tetra[ktet[1]].v[MMG5_idir[iface[1]][1]];
1375 ic1 = mesh->tetra[ktet[1]].v[MMG5_idir[iface[1]][2]];
1376
1377 assert ( ( (ia0 == ia1) && ((ib0 == ic1) && (ic0 == ib1 )) ) ||
1378 ( (ia0 == ib1) && ((ib0 == ia1) && (ic0 == ic1 )) ) ||
1379 ( (ia0 == ic1) && ((ib0 == ib1) && (ic0 == ia1 )) ) );
1380#endif
1381 }
1382
1383 return 1;
1384}
1385
1386
1388 double lmax, MMG5_int *eltab,int8_t metRidTyp) {
1389 MMG5_pTetra pt;
1390 MMG5_Hash hash;
1391 double len;
1392 MMG5_int k,np,nq;
1393 int8_t ia,i0,i1,ier;
1394
1395 /* Hash all edges in the mesh */
1396 if ( !MMG5_hashNew(mesh,&hash,mesh->np,7*mesh->np) ) return 0;
1397
1398 for(k=1; k<=mesh->ne; k++) {
1399 pt = &mesh->tetra[k];
1400 if ( !MG_EOK(pt) ) continue;
1401
1402 for(ia=0; ia<6; ia++) {
1403 i0 = MMG5_iare[ia][0];
1404 i1 = MMG5_iare[ia][1];
1405 np = pt->v[i0];
1406 nq = pt->v[i1];
1407
1408 if(!MMG5_hashEdge(mesh,&hash,np,nq,0)){
1409 fprintf(stderr,"\n ## Error: %s: function MMG5_hashEdge return 0\n",
1410 __func__);
1411 return 0;
1412 }
1413 }
1414 }
1415
1416 /* Pop edges from hash table, and analyze their length */
1417 for(k=1; k<=mesh->ne; k++) {
1418 pt = &mesh->tetra[k];
1419 if ( !MG_EOK(pt) ) continue;
1420
1421 for(ia=0; ia<6; ia++) {
1422 i0 = MMG5_iare[ia][0];
1423 i1 = MMG5_iare[ia][1];
1424 np = pt->v[i0];
1425 nq = pt->v[i1];
1426
1427 /* Remove edge from hash ; ier = 1 if edge has been found */
1428 ier = MMG5_hashPop(&hash,np,nq);
1429 if( ier ) {
1430 if ( (!metRidTyp) && met->m && met->size>1 ) {
1431 // Warning: for aniso metrics, we may erroneously approximate the
1432 // length of a curve boundary edge by the length of the straight edge
1433 // if the "MG_BDY" tag is missing along the edge.
1434 len = MMG5_lenedg(mesh,met,ia,pt);
1435 }
1436 else {
1437 // Warning: we may erroneously approximate the length of a curve
1438 // boundary edge by the length of the straight edge if the "MG_BDY"
1439 // tag is missing along the edge.
1440 len = MMG5_lenedg33_ani(mesh,met,ia,pt);
1441 }
1442
1443 if( (len < lmin) || (len > lmax) ) {
1444 eltab[k] = 1;
1445 break;
1446 }
1447 }
1448 }
1449 }
1450 MMG5_DEL_MEM(mesh,hash.item);
1451 return 1;
1452}
1453
1465static inline
1467 MMG5_pTetra pt;
1468 MMG5_int k;
1469 int i,ier;
1470
1471 assert ( mesh->info.optim );
1472
1473 /* Detect the point not used by the mesh */
1474 ++mesh->base;
1475
1476#ifndef NDEBUG
1477 for (k=1; k<=mesh->np; k++) {
1478 assert ( mesh->point[k].flag < mesh->base );
1479 }
1480#endif
1481
1482 /* For mmg3d, detect points used by tetra */
1483 for (k=1; k<=mesh->ne; k++) {
1484 pt = &mesh->tetra[k];
1485 if ( !MG_EOK(pt) ) continue;
1486
1487 for (i=0; i<4; i++) {
1488 mesh->point[pt->v[i]].flag = mesh->base;
1489 }
1490 }
1491
1492 /* Compute hmin/hmax on unflagged points and truncate the metric */
1493 if ( !ani ) {
1495 }
1496 else {
1497 MMG5_solTruncature_ani = MMG5_3dSolTruncature_ani;
1499 }
1500
1501 return ier;
1502}
1503
1515 MMG5_pTetra pt;
1516 MMG5_pPoint p1,p2;
1517 double ux,uy,uz,dd;
1518 MMG5_int k,ipa,ipb;
1519 int i,ia,ib,type;
1520 // we guess that we have less than INT32_MAX edges passing through a point
1521 int *mark;
1522
1523 MMG5_SAFE_CALLOC(mark,mesh->np+1,int,return 0);
1524
1525 /* Memory alloc */
1526 if ( met->size!=1 ) {
1527 fprintf(stderr,"\n ## Error: %s: unexpected size of metric: %d.\n",
1528 __func__,met->size);
1529 return 0;
1530 }
1531
1532 type=1;
1533 if ( !MMG3D_Set_solSize(mesh,met,MMG5_Vertex,mesh->np,type) )
1534 return 0;
1535
1536 /* Travel the triangles edges and add the edge contribution to edges
1537 * extermities */
1538 for (k=1; k<=mesh->ne; k++) {
1539 pt = &mesh->tetra[k];
1540 if ( !MG_EOK(pt) ) continue;
1541
1542 for (i=0; i<6; i++) {
1543 ia = MMG5_iare[i][0];
1544 ib = MMG5_iare[i][1];
1545 ipa = pt->v[ia];
1546 ipb = pt->v[ib];
1547 p1 = &mesh->point[ipa];
1548 p2 = &mesh->point[ipb];
1549
1550 ux = p1->c[0] - p2->c[0];
1551 uy = p1->c[1] - p2->c[1];
1552 uz = p1->c[2] - p2->c[2];
1553 dd = sqrt(ux*ux + uy*uy + uz*uz);
1554
1555 met->m[ipa] += dd;
1556 mark[ipa]++;
1557 met->m[ipb] += dd;
1558 mark[ipb]++;
1559 }
1560 }
1561
1562 /* vertex size */
1563 for (k=1; k<=mesh->np; k++) {
1564 if ( !mark[k] ) {
1565 continue;
1566 }
1567 else
1568 met->m[k] = met->m[k] / (double)mark[k];
1569 }
1570
1571 MMG5_SAFE_FREE(mark);
1572
1573 /* Metric truncation */
1575
1576 return 1;
1577}
1578
1590 MMG5_pTetra pt;
1591 MMG5_pPoint p1,p2;
1592 double u[3],dd,tensordot[6];
1593 MMG5_int k,iadr,ipa,ipb;
1594 int i,j,type;
1595 // we guess that we have less than INT32_MAX edges passing through a point
1596 int *mark;
1597
1598 MMG5_SAFE_CALLOC(mark,mesh->np+1,int,return 0);
1599
1600 /* Memory alloc */
1601 if ( met->size!=6 ) {
1602 fprintf(stderr,"\n ## Error: %s: unexpected size of metric: %d.\n",
1603 __func__,met->size);
1604 return 0;
1605 }
1606
1607 type = 3;
1608 if ( !MMG3D_Set_solSize(mesh,met,MMG5_Vertex,mesh->np,type) )
1609 return 0;
1610
1611 /* Travel the tetra edges and add the edge contribution to edges
1612 * extermities */
1613 for (k=1; k<=mesh->ne; k++) {
1614 pt = &mesh->tetra[k];
1615 if ( !MG_EOK(pt) ) continue;
1616
1617 for (i=0; i<6; i++) {
1618 ipa = pt->v[MMG5_iare[i][0]];
1619 ipb = pt->v[MMG5_iare[i][1]];
1620 p1 = &mesh->point[ipa];
1621 p2 = &mesh->point[ipb];
1622
1623 u[0] = p1->c[0] - p2->c[0];
1624 u[1] = p1->c[1] - p2->c[1];
1625 u[2] = p1->c[2] - p2->c[2];
1626
1627 tensordot[0] = u[0]*u[0];
1628 tensordot[1] = u[0]*u[1];
1629 tensordot[2] = u[0]*u[2];
1630 tensordot[3] = u[1]*u[1];
1631 tensordot[4] = u[1]*u[2];
1632 tensordot[5] = u[2]*u[2];
1633
1634 iadr = 6*ipa;
1635 for ( j=0; j<6; ++j ) {
1636 met->m[iadr+j] += tensordot[j];
1637 }
1638 mark[ipa]++;
1639
1640 iadr = 6*ipb;
1641 for ( j=0; j<6; ++j ) {
1642 met->m[iadr+j] += tensordot[j];
1643 }
1644 mark[ipb]++;
1645 }
1646 }
1647
1648 for (k=1; k<=mesh->np; k++) {
1649 if ( !mark[k] ) {
1650 continue;
1651 }
1652
1653 /* Metric = nedges/dim * inv (sum(tensor_dot(edges,edges))).
1654 * sum(tensor_dot) is stored in sol->m so reuse tensordot to
1655 * compute M. */
1656 iadr = 6*k;
1657 if ( !MMG5_invmat(met->m+iadr,tensordot) ) {
1658 /* Non invertible matrix: impose FLT_MIN, it will be truncated by hmax
1659 * later */
1660 fprintf(stdout, " ## Warning: %s: %d: non invertible matrix."
1661 " Impose hmax size at point\n",__func__,__LINE__);
1662 met->m[iadr+0] = FLT_MIN;
1663 met->m[iadr+1] = 0;
1664 met->m[iadr+2] = 0;
1665 met->m[iadr+3] = FLT_MIN;
1666 met->m[iadr+4] = 0;
1667 met->m[iadr+5] = FLT_MIN;
1668 continue;
1669 }
1670
1671 dd = (double)mark[k]/3.;
1672
1673 for ( j=0; j<6; ++j ) {
1674 met->m[iadr+j] = dd*tensordot[j];
1675 }
1676
1677#ifndef NDEBUG
1678 /* Check metric */
1679 double lambda[3],vp[3][3];
1680 if (!MMG5_eigenv3d(1,met->m+iadr,lambda,vp) ) {
1681 fprintf(stdout, " ## Warning: %s: %d: non diagonalizable metric.",
1682 __func__,__LINE__);
1683 }
1684
1685 assert ( lambda[0] > 0. && lambda[1] > 0. && lambda[2] > 0.
1686 && "Negative eigenvalue" );
1687 assert ( isfinite(lambda[0]) && isfinite(lambda[1]) && isfinite(lambda[2])
1688 && "Infinite eigenvalue" );
1689#endif
1690 }
1691
1692 MMG5_SAFE_FREE(mark);
1693
1695
1696 return 1;
1697}
1698
1700 double hsiz;
1701 int type;
1702
1703 /* Set solution size */
1704 if ( mesh->info.ani ) {
1705 met->size = 6;
1706 type = 3;
1707 }
1708 else {
1709 met->size = 1;
1710 type = 1;
1711 }
1712
1713 /* Memory alloc */
1714 if ( !MMG3D_Set_solSize(mesh,met,MMG5_Vertex,mesh->np,type) )
1715 return 0;
1716
1717 if ( !MMG5_Compute_constantSize(mesh,met,&hsiz) )
1718 return 0;
1719
1720 mesh->info.hsiz = hsiz;
1721
1722 MMG5_Set_constantSize(mesh,met,hsiz);
1723
1724 return 1;
1725}
1726
1728 MMG5_int k;
1729 double tmp;
1730
1731 if ( met->size!=6 ) { return 1; }
1732
1733 for ( k=1; k<=met->np; ++k ) {
1734
1735 tmp = met->m[6*k+2];
1736
1737 met->m[6*k+2] = met->m[6*k+3];
1738 met->m[6*k+3] = tmp;
1739
1740 }
1741 return 1;
1742}
1743
1744int MMG3D_Compute_eigenv(double m[6],double lambda[3],double vp[3][3]) {
1745
1746 return MMG5_eigenv3d(1,m,lambda,vp);
1747
1748}
1749
1751
1752 /* sol */
1753 if ( !sol ) return;
1754
1755 if ( sol->m )
1756 MMG5_DEL_MEM(mesh,sol->m);
1757
1758 if ( sol->namein ) {
1760 }
1761
1762 if ( sol->nameout ) {
1764 }
1765
1766 memset ( sol, 0x0, sizeof(MMG5_Sol) );
1767
1768 /* Reset state to a scalar status */
1769 sol->dim = 3;
1770 sol->ver = 2;
1771 sol->size = 1;
1772 sol->type = 1;
1773
1774 return;
1775}
1776
1778 MMG5_int k,nref;
1779
1780 nref = 0;
1782 if ( mesh->tria ) {
1783
1784 MMG5_int nt = mesh->nt;
1785 k = 1;
1786 do {
1787 MMG5_pTria ptt = &mesh->tria[k];
1788
1789 if ( !MG_EOK(ptt) ) {
1790 continue;
1791 }
1792
1793 if ( MMG5_abs(ptt->ref) == mesh->info.isoref ) {
1794 /* Current tria will be suppressed: search last non isosurf tria to fill
1795 * empty position */
1796 MMG5_pTria ptt1 = &mesh->tria[mesh->nt];
1797 assert( ptt1 );
1798
1799 while ( ((!MG_EOK(ptt1)) || (MMG5_abs(ptt1->ref) == mesh->info.isoref))
1800 && k < mesh->nt ) {
1801 --mesh->nt;
1802 ptt1 = &mesh->tria[mesh->nt];
1803
1804 }
1805 if ( ptt != ptt1 ) {
1806 /* We don't find any tria to keep after index k */
1807 memcpy(ptt,ptt1,sizeof(MMG5_Tria));
1808 --mesh->nt;
1809 }
1810 }
1811 /* Initially negative refs were used to mark isosurface: keep following
1812 * piece of code for retrocompatibility */
1813 if ( ptt->ref < 0 ) {
1814 ptt->ref = -ptt->ref;
1815 ++nref;
1816 }
1817 }
1818 while ( ++k < mesh->nt );
1819
1820 /* At the end of the loop, either k==mesh->nt, either k==mesh->nt+1 (because
1821 * tria at idx mesh->nt was iso or unused and element mesh->nt+1 has been
1822 * copied into k) */
1823 assert ( (k==mesh->nt) || (k==mesh->nt+1) );
1824
1825 /* Check if last element is iso */
1826 MMG5_pTria ptt = &mesh->tria[mesh->nt];
1827 if ( (!MG_EOK(ptt)) || (MMG5_abs(ptt->ref) == mesh->info.isoref) ) {
1828 --mesh->nt;
1829 }
1830
1831 if ( mesh->info.imprim > 4 ) {
1832 fprintf(stdout," Deleted iso triangles: %" MMG5_PRId "\n",nt-mesh->nt);
1833 }
1834
1835 if( !mesh->nt ) {
1837 }
1838 else if ( mesh->nt < nt ) {
1839 MMG5_ADD_MEM(mesh,(mesh->nt-nt)*sizeof(MMG5_Tria),"triangles",
1840 fprintf(stderr," Exit program.\n");
1841 return 0);
1843 "triangles",return 0);
1844 }
1845 }
1846
1848 return MMG5_Clean_isoEdges(mesh);
1849}
int MMG5_Compute_constantSize(MMG5_pMesh mesh, MMG5_pSol met, double *hsiz)
void MMG5_Set_constantSize(MMG5_pMesh mesh, MMG5_pSol met, double hsiz)
char * MMG5_Get_filenameExt(char *filename)
int MMG3D_Set_inputMeshName(MMG5_pMesh mesh, const char *meshin)
Set the name of input mesh.
int MMG3D_Set_inputParamName(MMG5_pMesh mesh, const char *fparamin)
Set the name of the input parameter file.
int MMG3D_Set_dparameter(MMG5_pMesh mesh, MMG5_pSol sol, int dparam, double val)
set a real-valued parameter of the remesher
int MMG3D_Set_localParameter(MMG5_pMesh mesh, MMG5_pSol sol, int typ, MMG5_int ref, double hmin, double hmax, double hausd)
set a local parameter
int MMG3D_Set_lsBaseReference(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int br)
Set a new level-set base reference.
int MMG3D_Set_inputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solin)
Set the name of input solution file.
int MMG3D_Set_multiMat(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int ref, int split, MMG5_int rmin, MMG5_int rplus)
Set the reference mapping for the elements of reference ref in level-set discretization mode.
int MMG3D_Set_outputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solout)
Set the name of the output solution file.
int MMG3D_Set_outputMeshName(MMG5_pMesh mesh, const char *meshout)
Set the name of output mesh file.
int MMG3D_Set_iparameter(MMG5_pMesh mesh, MMG5_pSol sol, int iparam, MMG5_int val)
set an integer parameter of the remesher
int MMG3D_Set_solSize(MMG5_pMesh mesh, MMG5_pSol sol, int typEntity, MMG5_int np, int typSol)
Initialize a solution field.
int ier
tmp[*strlen0]
MMG5_pMesh MMG5_pSol * sol
if(!ier) exit(EXIT_FAILURE)
MMG5_pMesh * mesh
int MMG3D_PROctreein_iso(MMG5_pMesh mesh, MMG5_pSol sol, MMG3D_pPROctree PROctree, MMG5_int ip, double lmax)
Definition: PRoctree_3d.c:1143
int MMG3D_PROctreein_ani(MMG5_pMesh mesh, MMG5_pSol sol, MMG3D_pPROctree PROctree, MMG5_int ip, double lmax)
Definition: PRoctree_3d.c:1225
int MMG5_setadj(MMG5_pMesh mesh)
Definition: analys_3d.c:108
int MMG5_movbdynompt_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *listv, int ilistv, MMG5_int *lists, int ilists, int improve)
int MMG5_movbdyregpt_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *listv, int ilistv, MMG5_int *lists, int ilists, int improveSurf, int improveVol)
int MMG5_movbdyrefpt_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *listv, int ilistv, MMG5_int *lists, int ilists, int improve)
int MMG5_movintpt_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *list, int ilist, int improve)
Definition: anisomovpt_3d.c:58
int MMG5_movbdyridpt_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG3D_pPROctree PROctree, int64_t *listv, int ilistv, MMG5_int *lists, int ilists, int improve)
int MMG5_compute_meanMetricAtMarkedPoints_ani(MMG5_pMesh mesh, MMG5_pSol met)
Definition: anisosiz.c:2205
int MMG3D_gradsiz_ani(MMG5_pMesh mesh, MMG5_pSol met)
Definition: anisosiz_3d.c:2073
int MMG3D_defsiz_ani(MMG5_pMesh mesh, MMG5_pSol met)
Definition: anisosiz_3d.c:1310
int MMG3D_gradsizreq_ani(MMG5_pMesh mesh, MMG5_pSol met)
Definition: anisosiz_3d.c:2210
MMG5_Info info
void tminit(mytime *t, int maxtim)
Initialize mytime object.
Definition: chrono.c:120
void printim(double elps, char *stim)
Print real time.
Definition: chrono.c:160
void chrono(int cmode, mytime *ptt)
Function to measure time.
Definition: chrono.c:49
#define TIMEMAX
int MMG5_cavity_iso(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int iel, int ip, int64_t *list, int lon, double volmin)
Definition: delaunay_3d.c:741
int MMG5_cavity_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int iel, int ip, int64_t *list, int lon, double volmin)
Definition: delaunay_3d.c:581
int MMG5_eigenv3d(int symmat, double *mat, double lambda[3], double v[3][3])
Find eigenvalues and vectors of a 3x3 matrix.
Definition: eigenv.c:385
int MMG5_hashNew(MMG5_pMesh mesh, MMG5_Hash *hash, MMG5_int hsiz, MMG5_int hmax)
Definition: hash.c:532
int MMG5_hashEdge(MMG5_pMesh mesh, MMG5_Hash *hash, MMG5_int a, MMG5_int b, MMG5_int k)
Definition: hash.c:346
int MMG3D_hashTetra(MMG5_pMesh mesh, int pack)
Create array of adjacency.
Definition: hash_3d.c:122
int MMG3D_hashPrism(MMG5_pMesh mesh)
Definition: hash_3d.c:240
int MMG5_chkBdryTria(MMG5_pMesh mesh)
Definition: hash_3d.c:1559
int MMG5_bdrySet(MMG5_pMesh mesh)
Definition: hash_3d.c:1958
int MMG3D_hashTria(MMG5_pMesh mesh, MMG5_Hash *hash)
Definition: hash_3d.c:838
int MMG5_hashPop(MMG5_Hash *hash, MMG5_int a, MMG5_int b)
Definition: hash_3d.c:850
int MMG5_bdryPerm(MMG5_pMesh mesh)
Definition: hash_3d.c:2385
static double MMG5_caltet_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTetra pt)
static double MMG5_caltet_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTetra pt)
static double MMG5_lenedg_iso(MMG5_pMesh mesh, MMG5_pSol met, int ia, MMG5_pTetra pt)
static double MMG5_lenedg33_ani(MMG5_pMesh mesh, MMG5_pSol met, int ia, MMG5_pTetra pt)
static double MMG3D_caltetLES_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTetra pt)
static double MMG5_lenedg_ani(MMG5_pMesh mesh, MMG5_pSol met, int ia, MMG5_pTetra pt)
static double MMG5_lenedgCoor_ani(double *ca, double *cb, double *sa, double *sb)
inlined Functions
static double MMG5_lenSurfEdg_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int ip1, MMG5_int ip2, int8_t isedg)
static double MMG5_lenSurfEdg_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int np0, MMG5_int np1, int8_t isedg)
int MMG5_interp4bar_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int k, MMG5_int ip, double cb[4])
Definition: intmet_3d.c:450
int MMG5_interp4bar_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int k, MMG5_int ip, double cb[4])
Definition: intmet_3d.c:380
int MMG5_intmet_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int k, int8_t i, MMG5_int ip, double s)
Definition: intmet_3d.c:174
int MMG5_intmet_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int k, int8_t i, MMG5_int ip, double s)
Definition: intmet_3d.c:57
int MMG5_compute_meanMetricAtMarkedPoints_iso(MMG5_pMesh mesh, MMG5_pSol met)
Definition: isosiz.c:167
double MMG5_lenedgCoor_iso(double *ca, double *cb, double *ma, double *mb)
Compute edge length from edge's coordinates.
Definition: isosiz_3d.c:61
int MMG3D_gradsizreq_iso(MMG5_pMesh mesh, MMG5_pSol met)
Definition: isosiz_3d.c:1168
int MMG3D_gradsiz_iso(MMG5_pMesh mesh, MMG5_pSol met)
Definition: isosiz_3d.c:1083
int MMG3D_defsiz_iso(MMG5_pMesh mesh, MMG5_pSol met)
Definition: isosiz_3d.c:662
void MMG3D_Free_topoTables(MMG5_pMesh mesh)
Definition: libmmg3d.c:67
void MMG3D_Set_commonFunc(void)
Definition: libmmg3d.c:1745
API headers and documentation for the mmg3d library, for volumetric meshes in 3D.
LIBMMG3D_EXPORT double(* MMG3D_lenedgCoor)(double *ca, double *cb, double *sa, double *sb)
Compute the length of an edge according to the size prescription.
Definition: mmg3dexterns.c:10
LIBMMG3D_EXPORT int(* MMG3D_doSol)(MMG5_pMesh mesh, MMG5_pSol met)
Compute isotropic size map according to the mean of the length of the edges passing through a vertex.
Definition: mmg3dexterns.c:11
@ MMG3D_DPARAM_hmin
Definition: libmmg3d.h:171
@ MMG3D_IPARAM_debug
Definition: libmmg3d.h:146
@ MMG3D_IPARAM_numberOfLocalParam
Definition: libmmg3d.h:161
@ MMG3D_IPARAM_isoref
Definition: libmmg3d.h:169
@ MMG3D_IPARAM_noswap
Definition: libmmg3d.h:156
@ MMG3D_IPARAM_opnbdy
Definition: libmmg3d.h:151
@ MMG3D_DPARAM_angleDetection
Definition: libmmg3d.h:170
@ MMG3D_DPARAM_hsiz
Definition: libmmg3d.h:173
@ MMG3D_IPARAM_isosurf
Definition: libmmg3d.h:149
@ MMG3D_DPARAM_hausd
Definition: libmmg3d.h:174
@ MMG3D_IPARAM_nosurf
Definition: libmmg3d.h:158
@ MMG3D_IPARAM_nomove
Definition: libmmg3d.h:157
@ MMG3D_IPARAM_renum
Definition: libmmg3d.h:165
@ MMG3D_DPARAM_hgrad
Definition: libmmg3d.h:175
@ MMG3D_IPARAM_nosizreq
Definition: libmmg3d.h:168
@ MMG3D_DPARAM_rmc
Definition: libmmg3d.h:179
@ MMG3D_IPARAM_angle
Definition: libmmg3d.h:147
@ MMG3D_IPARAM_noinsert
Definition: libmmg3d.h:155
@ MMG3D_DPARAM_ls
Definition: libmmg3d.h:177
@ MMG3D_DPARAM_hgradreq
Definition: libmmg3d.h:176
@ MMG3D_IPARAM_xreg
Definition: libmmg3d.h:160
@ MMG3D_DPARAM_hmax
Definition: libmmg3d.h:172
@ MMG3D_IPARAM_numsubdomain
Definition: libmmg3d.h:164
@ MMG3D_IPARAM_lag
Definition: libmmg3d.h:152
@ MMG3D_IPARAM_nreg
Definition: libmmg3d.h:159
@ MMG3D_IPARAM_numberOfMat
Definition: libmmg3d.h:163
@ MMG3D_IPARAM_verbose
Definition: libmmg3d.h:144
@ MMG3D_DPARAM_xreg
Definition: libmmg3d.h:178
@ MMG3D_IPARAM_numberOfLSBaseReferences
Definition: libmmg3d.h:162
@ MMG3D_IPARAM_optimLES
Definition: libmmg3d.h:154
@ MMG3D_IPARAM_optim
Definition: libmmg3d.h:153
@ MMG3D_IPARAM_iso
Definition: libmmg3d.h:148
@ MMG3D_IPARAM_octree
Definition: libmmg3d.h:167
@ MMG3D_IPARAM_nofem
Definition: libmmg3d.h:150
@ MMG3D_IPARAM_mem
Definition: libmmg3d.h:145
#define MMG3D_ALPHAD
int MMG5_movbdyregpt_iso(MMG5_pMesh, MMG5_pSol, MMG3D_pPROctree, int64_t *, int, MMG5_int *, int, int, int)
Definition: movpt_3d.c:626
void MMG5_freeXPrisms(MMG5_pMesh mesh)
Definition: zaldy_3d.c:378
void MMG5_freeXTets(MMG5_pMesh mesh)
Definition: zaldy_3d.c:359
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
int MMG5_movbdyrefpt_iso(MMG5_pMesh, MMG5_pSol, MMG3D_pPROctree, int64_t *, int, MMG5_int *, int, int)
Definition: movpt_3d.c:1446
int MMG5_movintpt_iso(MMG5_pMesh, MMG5_pSol, MMG3D_pPROctree, int64_t *, int, int)
Definition: movpt_3d.c:58
static const uint8_t MMG5_idir_pr[5][4]
idir[i]: vertices of face i for a prism
int MMG5_movbdyridpt_iso(MMG5_pMesh, MMG5_pSol, MMG3D_pPROctree, int64_t *, int, MMG5_int *, int, int)
Definition: movpt_3d.c:1614
int MMG3D_memOption(MMG5_pMesh mesh)
Definition: zaldy_3d.c:271
static const uint8_t MMG5_idir[4][3]
idir[i]: vertices of face opposite to vertex i
int MMG5_movbdynompt_iso(MMG5_pMesh, MMG5_pSol, MMG3D_pPROctree, int64_t *, int, MMG5_int *, int, int)
Definition: movpt_3d.c:1470
static void MMG5_warnOrientation(MMG5_pMesh mesh)
#define MMG_ARGV_APPEND(fromV, toV, fromC, toC, on_failure)
void MMG3D_searchqua(MMG5_pMesh mesh, MMG5_pSol met, double critmin, MMG5_int *eltab, int8_t metRidTyp)
List bad elements.
int MMG3D_searchlen(MMG5_pMesh mesh, MMG5_pSol met, double lmin, double lmax, MMG5_int *eltab, int8_t metRidTyp)
List edges that are too short or too long.
void MMG3D_destockOptions(MMG5_pMesh mesh, MMG5_Info *info)
Recover the info structure stored in the mesh structure.
int MMG3D_usage(char *prog)
Print help for mmg3d options.
static int MMG3D_solTruncatureForOptim(MMG5_pMesh mesh, MMG5_pSol met, int ani)
void MMG3D_setfunc(MMG5_pMesh mesh, MMG5_pSol met)
Set function pointers for caltet, lenedg, lenedgCoor defsiz, gradsiz... depending if the metric that ...
int MMG3D_mmg3dcheck(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, double critmin, double lmin, double lmax, MMG5_int *eltab, int8_t metRidTyp)
Search invalid elements (in term of quality or edge length) in a mesh.
int MMG3D_Get_nonBdyTriangle(MMG5_pMesh mesh, MMG5_int *v0, MMG5_int *v1, MMG5_int *v2, MMG5_int *ref, MMG5_int idx)
Get vertices and reference of a non-boundary triangle.
int MMG3D_parsop(MMG5_pMesh mesh, MMG5_pSol met)
Read a file containing Local parameters (.mmg3d extension)
void MMG5_argv_cleanup(char **mmgArgv, int mmgArgc)
int MMG3D_Compute_eigenv(double m[6], double lambda[3], double vp[3][3])
Compute the real eigenvalues and eigenvectors of a symmetric matrix.
int MMG3D_Get_tetsFromTria(MMG5_pMesh mesh, MMG5_int ktri, MMG5_int ktet[2], int iface[2])
Get two tetrahedra given a triangle and face indices.
void MMG3D_Free_solutions(MMG5_pMesh mesh, MMG5_pSol sol)
Free the solution structure of a given mesh.
int MMG3D_Get_numberOfNonBdyTriangles(MMG5_pMesh mesh, MMG5_int *nb_tria)
Get the number of non-boundary triangles.
int MMG3D_stockOptions(MMG5_pMesh mesh, MMG5_Info *info)
Store the info structure in the mesh structure.
int MMG3D_freeLocalPar(MMG5_pMesh mesh)
int MMG3D_doSol_iso(MMG5_pMesh mesh, MMG5_pSol met)
int MMG3D_storeknownar(int argc, char *argv[], MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, int *mmgArgc, char *mmgArgv[])
int MMG3D_Set_constantSize(MMG5_pMesh mesh, MMG5_pSol met)
Compute a constant size map according to the hsiz, hmin and hmax parameters.
int MMG3D_defaultValues(MMG5_pMesh mesh)
Print the default parameters values.
int MMG3D_Clean_isoSurf(MMG5_pMesh mesh)
Clean data (triangles and edges) linked to isosurface.
int MMG3D_Get_tetFromTria(MMG5_pMesh mesh, MMG5_int ktri, MMG5_int *ktet, int *iface)
Get a tetrahedron given one of its triangles and the index by which it refers to this triangle (DEPRE...
int MMG3D_parsar(int argc, char *argv[], MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol)
Store command-line arguments.
int MMG3D_Get_adjaTet(MMG5_pMesh mesh, MMG5_int kel, MMG5_int listet[4])
Return adjacent elements of a tetrahedron.
int MMG3D_switch_metricStorage(MMG5_pMesh mesh, MMG5_pSol met)
Swap the m22 and m23 values of the metric.
int MMG3D_doSol_ani(MMG5_pMesh mesh, MMG5_pSol met)
LIBMMG_CORE_EXPORT int MMG5_Clean_isoEdges(MMG5_pMesh mesh)
Definition: libtools.c:368
LIBMMG_CORE_EXPORT int MMG5_scaleMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol ls)
Definition: scalem.c:649
#define MMG5_MMAT_Split
Definition: libmmgtypes.h:211
#define MMG5_STRONGFAILURE
Definition: libmmgtypes.h:65
#define MMG5_LOWFAILURE
Definition: libmmgtypes.h:57
#define MMG5_SUCCESS
Definition: libmmgtypes.h:49
@ MMG5_Tensor
Definition: libmmgtypes.h:221
#define MMG5_MMAT_NoSplit
Definition: libmmgtypes.h:203
@ MMG5_Vertex
Definition: libmmgtypes.h:230
@ MMG5_Tetrahedron
Definition: libmmgtypes.h:233
@ MMG5_Triangle
Definition: libmmgtypes.h:232
void MMG5_advancedUsage(void)
Definition: libtools.c:305
void MMG5_mmgUsage(char *prog)
Definition: libtools.c:212
void MMG5_paramUsage1(void)
Definition: libtools.c:242
void MMG5_2d3dUsage(void)
Definition: libtools.c:294
void MMG5_mmgDefaultValues(MMG5_pMesh mesh)
Definition: libtools.c:70
void MMG5_lagUsage(void)
Definition: libtools.c:278
void MMG5_paramUsage2(void)
Definition: libtools.c:261
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
int MMG5_3dSolTruncature_ani(MMG5_pMesh mesh, MMG5_pSol met)
Definition: scalem.c:453
#define _LIBMMG5_RETURN(mesh, sol, met, val)
double MMG5_caltri_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTria ptt)
Definition: quality.c:115
#define MG_EOK(pt)
static void MMG5_excfun(int sigid)
#define MMG5_ADD_MEM(mesh, size, message, law)
#define MMG5_UNSET
static void MMG5_warnScotch(MMG5_pMesh mesh)
int MMG5_solTruncature_iso(MMG5_pMesh mesh, MMG5_pSol met)
Definition: scalem.c:299
#define MG_BDY
#define MMG5_SAFE_MALLOC(ptr, size, type, law)
#define MG_STR
#define MMG5_FILESTR_LGTH
double MMG5_caltri_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pTria ptt)
Definition: quality.c:199
#define MMG5_SAFE_FREE(ptr)
#define MMG_FSCANF(stream, format,...)
#define MMG5_DEL_MEM(mesh, ptr)
int MMG5_invmat(double *m, double *mi)
Definition: tools.c:562
#define MMG5_SAFE_RECALLOC(ptr, prevSize, newSize, type, message, law)
Identic as MMG5_HGeom but use MMG5_hedge to store edges instead of MMG5_hgeom (memory economy).
Definition: libmmgtypes.h:603
MMG5_hedge * item
Definition: libmmgtypes.h:605
Structure to store input parameters of the job.
Definition: libmmgtypes.h:523
int8_t iso
Definition: libmmgtypes.h:541
double hsiz
Definition: libmmgtypes.h:525
int8_t isosurf
Definition: libmmgtypes.h:542
uint8_t ani
Definition: libmmgtypes.h:553
MMG5_int isoref
Definition: libmmgtypes.h:528
int PROctree
Definition: libmmgtypes.h:534
uint8_t optimLES
Definition: libmmgtypes.h:553
int8_t lag
Definition: libmmgtypes.h:547
char * fparam
Definition: libmmgtypes.h:555
MMG5_pPar par
Definition: libmmgtypes.h:524
uint8_t optim
Definition: libmmgtypes.h:553
MMG mesh structure.
Definition: libmmgtypes.h:613
MMG5_int ntmax
Definition: libmmgtypes.h:620
MMG5_Info info
Definition: libmmgtypes.h:659
MMG5_int * adjapr
Definition: libmmgtypes.h:640
MMG5_int ne
Definition: libmmgtypes.h:620
MMG5_pPoint point
Definition: libmmgtypes.h:649
char * nameout
Definition: libmmgtypes.h:661
MMG5_int mark
Definition: libmmgtypes.h:626
MMG5_int * adja
Definition: libmmgtypes.h:632
MMG5_pPrism prism
Definition: libmmgtypes.h:653
MMG5_int nemax
Definition: libmmgtypes.h:620
MMG5_int npmax
Definition: libmmgtypes.h:620
MMG5_int base
Definition: libmmgtypes.h:624
MMG5_pTetra tetra
Definition: libmmgtypes.h:651
MMG5_pxPrism xprism
Definition: libmmgtypes.h:654
MMG5_int nt
Definition: libmmgtypes.h:620
MMG5_pTria tria
Definition: libmmgtypes.h:655
MMG5_int np
Definition: libmmgtypes.h:620
MMG5_pxTetra xtetra
Definition: libmmgtypes.h:652
MMG5_int nprism
Definition: libmmgtypes.h:621
char * namein
Definition: libmmgtypes.h:660
Structure to store vertices of an MMG mesh.
Definition: libmmgtypes.h:276
double c[3]
Definition: libmmgtypes.h:277
MMG5_int flag
Definition: libmmgtypes.h:288
Structure to store prsim of a MMG mesh.
Definition: libmmgtypes.h:469
MMG5_int v[6]
Definition: libmmgtypes.h:470
MMG5_int ref
Definition: libmmgtypes.h:471
MMG5_int xpr
Definition: libmmgtypes.h:474
char * nameout
Definition: libmmgtypes.h:683
char * namein
Definition: libmmgtypes.h:682
double * m
Definition: libmmgtypes.h:680
MMG5_int np
Definition: libmmgtypes.h:674
Structure to store tetrahedra of an MMG mesh.
Definition: libmmgtypes.h:407
MMG5_int v[4]
Definition: libmmgtypes.h:409
MMG5_int xt
Definition: libmmgtypes.h:413
MMG5_int ref
Definition: libmmgtypes.h:410
Structure to store triangles of a MMG mesh.
Definition: libmmgtypes.h:338
MMG5_int ref
Definition: libmmgtypes.h:341
MMG5_int v[3]
Definition: libmmgtypes.h:340
MMG5_int cc
Definition: libmmgtypes.h:343
uint16_t ftag[5]
Definition: libmmgtypes.h:491
MMG5_int ref[5]
Definition: libmmgtypes.h:485
MMG5_int ref[4]
Definition: libmmgtypes.h:426
uint16_t ftag[4]
Definition: libmmgtypes.h:430
Chrono object.