Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
libmmg2d_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 "libmmg2d.h"
34#include "libmmg2d_private.h"
36#include "mmgexterns_private.h"
37
39 if ( mesh->info.ani || (met && met->size==3 ) ) {
40 /* Force data consistency: if aniso metric is provided, met->size==3 and
41 * info.ani==0; with -A option, met->size==1 and info.ani==1 */
42 met->size = 3;
43 mesh->info.ani = 1;
44
45 /* Set pointers */
46 MMG2D_lencurv = MMG2D_lencurv_ani;
47 MMG5_compute_meanMetricAtMarkedPoints = MMG5_compute_meanMetricAtMarkedPoints_ani;
48 MMG2D_defsiz = MMG2D_defsiz_ani;
49 MMG2D_gradsiz = lissmet_ani;
50 MMG2D_gradsizreq = MMG5_gradsizreq_ani;
51 MMG2D_caltri = MMG2D_caltri_ani;
52 MMG2D_intmet = MMG2D_intmet_ani;
54 }
55 else {
56 MMG2D_lencurv = MMG2D_lencurv_iso;
57 MMG5_compute_meanMetricAtMarkedPoints = MMG5_compute_meanMetricAtMarkedPoints_iso;
58 MMG2D_defsiz = MMG2D_defsiz_iso;
59 MMG2D_gradsiz = MMG5_gradsiz_iso;
60 MMG2D_gradsizreq = MMG5_gradsizreq_iso;
61 MMG2D_caltri = MMG2D_caltri_iso;
62 MMG2D_intmet = MMG2D_intmet_iso;
64 }
65 return;
66}
67
68int MMG2D_usage(char *name) {
69
70 /* Common generic options, file options and mode options */
71 MMG5_mmgUsage(name);
72
73 /* Lagrangian option (only for mmg2d/3d) */
75
76 /* Common parameters (first section) */
78
79 /* Parameters shared by mmg2d and 3d only*/
81
82 /* Specific parameters */
83 fprintf(stdout,"-3dMedit val read and write for gmsh visu: output only if val=1, input and output if val=2, input if val=3\n");
84 fprintf(stdout,"\n");
85
86 fprintf(stdout,"-nofem do not force Mmg to create a finite element mesh \n");
87 fprintf(stdout,"-nosurf no surface modifications\n");
88
89 /* Common parameters (second section) */
91
92 /* Common options for advanced users */
94
95 fprintf(stdout,"\n\n");
96
97 return 1;
98}
99
100// In ls mode : metric must be provided using -met option (-sol or default is the ls).
101// In adp mode : -sol or -met or default allow to store the metric.
102int MMG2D_parsar(int argc,char *argv[],MMG5_pMesh mesh,MMG5_pSol met,MMG5_pSol sol) {
103 MMG5_pSol tmp = NULL;
104 double val;
105 int i;
106 char namein[MMG5_FILESTR_LGTH],*endptr;
107 int param;
108
109 /* First step: search if user want to see the default parameters values. */
110 for ( i=1; i< argc; ++i ) {
111 if ( !strcmp(argv[i],"-val") ) {
113 return 0;
114 }
115 else if ( ( !strcmp( argv[ i ],"-?" ) ) || ( !strcmp( argv[ i ],"-h" ) ) ) {
116 MMG2D_usage(argv[0]);
117 return 0;
118 }
119 }
120
121 /* Second step: read all other arguments. */
122 i = 1;
123 while ( i < argc ) {
124 if ( *argv[i] == '-' ) {
125 switch(argv[i][1]) {
126 case 'a':
127 if ( !strcmp(argv[i],"-ar") ) {
128 if ( i >= argc -1 ) {
129 fprintf(stderr,"\nMissing argument option %s\n",argv[i]);
130 return 0;
131 }
132 else {
133 val = strtof(argv[i+1],&endptr);
134 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
135 ++i;
137 return 0;
138 }
139 else {
140 /* argument is not a number */
141 fprintf(stderr,"\nMissing argument option %s\n",argv[i]);
142 return 0;
143 }
144 }
145 }
146 break;
147 case 'A': /* anisotropy */
149 return 0;
150 break;
151 case 'd':
152 if ( !strcmp(argv[i],"-default") ) {
153 mesh->mark=1;
154 } else { /* debug */
156 return 0;
157 }
158 break;
159 case 'f':
160 if ( !strcmp(argv[i],"-f") ) {
161 if ( ++i < argc && isascii(argv[i][0]) && argv[i][0]!='-' ) {
162 if ( !MMG2D_Set_inputParamName(mesh,argv[i]) )
163 return 0;
164 }
165 else {
166 fprintf(stderr,"\nMissing filename for %s\n",argv[i-1]);
167 MMG2D_usage(argv[0]);
168 return 0;
169 }
170 }
171 break;
172 case 'h':
173 param = MMG5_UNSET;
174 if ( i >= argc -1 ) {
175 fprintf(stderr,"\nMissing argument option %s\n",argv[i]);
176 return 0;
177 }
178 else {
179 if ( !strcmp(argv[i],"-hmin") ) {
180 param = MMG2D_DPARAM_hmin;
181 }
182 else if ( !strcmp(argv[i],"-hmax") ) {
183 param = MMG2D_DPARAM_hmax;
184 }
185 else if ( !strcmp(argv[i],"-hsiz") ) {
186 param = MMG2D_DPARAM_hsiz;
187 }
188 else if ( !strcmp(argv[i],"-hausd") ) {
189 param = MMG2D_DPARAM_hausd;
190 }
191 else if ( !strcmp(argv[i],"-hgradreq") ) {
192 param = MMG2D_DPARAM_hgradreq;
193 }
194 else if ( !strcmp(argv[i],"-hgrad") ) {
195 param = MMG2D_DPARAM_hgrad;
196 }
197 else {
198 /* Arg unknown by Mmg: arg starts with -h but is not known */
199 MMG2D_usage(argv[0]);
200 return 0;
201 }
202
203 assert ( param != MMG5_UNSET );
204
205 val = strtof(argv[i+1],&endptr);
206 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
207 ++i;
208 if ( !MMG2D_Set_dparameter(mesh,met,param,val) ){
209 return 0;
210 }
211 } else {
212 /* argument is not a number */
213 fprintf(stderr,"\nMissing argument option %s\n",argv[i]);
214 return 0;
215 }
216 }
217
218 break;
219 case 'i':
220 if ( !strcmp(argv[i],"-in") ) {
221 if ( ++i < argc && isascii(argv[i][0]) && argv[i][0]!='-') {
222 if ( !MMG2D_Set_inputMeshName(mesh, argv[i]) )
223 return 0;
224
226 return 0;
227 }else{
228 fprintf(stderr,"\nMissing filname for %s\n",argv[i]);
229 MMG2D_usage(argv[0]);
230 return 0;
231 }
232 }
233 else if ( !strcmp(argv[i],"-isoref") && ++i <= argc ) {
235 atoi(argv[i])) )
236 return 0;
237 }
238 else {
239 MMG2D_usage(argv[0]);
240 return 0;
241 }
242 break;
243 case 'l':
244 if ( !strcmp(argv[i],"-lag") ) {
245 if ( ++i < argc && isdigit(argv[i][0]) ) {
246 if ( !MMG2D_Set_iparameter(mesh,met,MMG2D_IPARAM_lag,atoi(argv[i])) )
247 return 0;
248 }
249 else {
250 fprintf(stderr,"\nMissing or unexpected argument option %s\n",argv[i-1]);
251 MMG2D_usage(argv[0]);
252 return 0;
253 }
254 }
255 else if ( !strcmp(argv[i],"-ls") ) {
257 return 0;
258
259 if ( i < argc -1 ) {
260 val = strtof(argv[i+1],&endptr);
261 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
262 ++i;
264 return 0;
265 }
266 }
267 }
268 else if ( !strcmp(argv[i],"-lssurf") ) {
270 return 0;
271
272 if ( i < argc -1 ) {
273 val = strtof(argv[i+1],&endptr);
274 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
275 ++i;
277 return 0;
278 }
279 }
280 }
281 break;
282 case 'm': /* memory */
283 if ( !strcmp(argv[i],"-met") ) {
284 if ( !met ) {
285 fprintf(stderr,"\nNo metric structure allocated for %s option\n",
286 argv[i]);
287 return 0;
288 }
289 if ( ++i < argc && isascii(argv[i][0]) && argv[i][0]!='-' ) {
290 if ( !MMG2D_Set_inputSolName(mesh,met,argv[i]) )
291 return 0;
292 }
293 else {
294 fprintf(stderr,"\nMissing filname for %s\n",argv[i-1]);
295 MMG2D_usage(argv[0]);
296 return 0;
297 }
298 }
299 else if (!strcmp(argv[i],"-m") ) {
300 if ( ++i < argc && isdigit(argv[i][0]) ) {
301 if ( !MMG2D_Set_iparameter(mesh,met,MMG2D_IPARAM_mem,atoi(argv[i])) )
302 return 0;
303 }
304 else {
305 fprintf(stderr,"\nMissing argument option %s\n",argv[i-1]);
306 MMG2D_usage(argv[0]);
307 return 0;
308 }
309 }
310 break;
311 case 'n':
312 if ( !strcmp(argv[i],"-nofem") ) {
314 return 0;
315 }
316 if ( !strcmp(argv[i],"-nreg") ) {
318 return 0;
319 }
320 else if ( !strcmp(argv[i],"-nr") ) {
322 return 0;
323 }
324 else if ( !strcmp(argv[i],"-nsd") ) {
325 if ( ++i < argc && isdigit(argv[i][0]) ) {
326 if ( !MMG2D_Set_iparameter(mesh,met,MMG2D_IPARAM_numsubdomain,atoi(argv[i])) )
327 return 0;
328 }
329 else {
330 fprintf(stderr,"\nMissing argument option %s\n",argv[i-1]);
331 MMG2D_usage(argv[0]);
332 return 0;
333 }
334 } else if ( !strcmp(argv[i],"-noswap") ) {
336 return 0;
337 }
338 else if( !strcmp(argv[i],"-noinsert") ) {
340 return 0;
341 }
342 else if( !strcmp(argv[i],"-nomove") ) {
344 return 0;
345 }
346 else if( !strcmp(argv[i],"-nosurf") ) {
348 return 0;
349 }
350 else if( !strcmp(argv[i],"-nosizreq") ) {
352 return 0;
353 }
354 }
355 break;
356 case 'o':
357 if ( (!strcmp(argv[i],"-out")) || (!strcmp(argv[i],"-o")) ) {
358 if ( ++i < argc && isascii(argv[i][0]) && argv[i][0]!='-') {
359 if ( !MMG2D_Set_outputMeshName(mesh,argv[i]) )
360 return 0;
361 }else{
362 fprintf(stderr,"\nMissing filname for %s\n",argv[i-1]);
363 MMG2D_usage(argv[0]);
364 return 0;
365 }
366 }
367 else if ( !strcmp(argv[i],"-opnbdy") ) {
369 return 0;
370 }
371 else if( !strcmp(argv[i],"-optim") ) {
373 return 0;
374 }
375 break;
376 case 'r':
377 if ( !strcmp(argv[i],"-rmc") ) {
379 return 0;
380 if ( i < argc -1 ) {
381 val = strtof(argv[i+1],&endptr);
382 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
383 ++i;
385 return 0;
386 }
387 }
388 }
389 break;
390 case 's':
391 if ( !strcmp(argv[i],"-sol") ) {
392 /* For retrocompatibility, store the metric if no sol structure available */
393 tmp = sol ? sol : met;
394
395 assert(tmp);
396 if ( ++i < argc && isascii(argv[i][0]) && argv[i][0]!='-' ) {
397 if ( !MMG2D_Set_inputSolName(mesh,tmp,argv[i]) )
398 return 0;
399 }
400 else {
401 fprintf(stderr,"\nMissing filname for %s\n",argv[i-1]);
402 MMG2D_usage(argv[0]);
403 return 0;
404 }
405 }
406 break;
407 case 'v':
408 if ( ++i < argc ) {
409 if ( argv[i][0] == '-' || isdigit(argv[i][0]) ) {
410 if ( !MMG2D_Set_iparameter(mesh,met,MMG2D_IPARAM_verbose,atoi(argv[i])) )
411 return 0;
412 }
413 else
414 i--;
415 }
416 else {
417 fprintf(stderr,"\nMissing argument option for %s\n",argv[i-1]);
418 MMG2D_usage(argv[0]);
419 return 0;
420 }
421 break;
422 case 'x':
423 if ( !strcmp(argv[i],"-xreg") ) {
425 return 0;
426 if ( i < argc -1 ) {
427 val = strtof(argv[i+1],&endptr);
428 if ( endptr == &(argv[i+1][strlen(argv[i+1])]) ) {
429 ++i;
431 return 0;
432 }
433 }
434 }
435 break;
436 case '3':
437 if(!strcmp(argv[i],"-3dMedit") ) {
438 if ( ++i < argc && isdigit(argv[i][0]) ) {
439 if ( !MMG2D_Set_iparameter(mesh,met,MMG2D_IPARAM_3dMedit,atoi(argv[i])) )
440 return 0;
441 }
442 else {
443 fprintf(stderr,"\nMissing argument option for %s\n",argv[i-1]);
444 MMG2D_usage(argv[0]);
445 return 0;
446 }
447 }
448 break;
449 default:
450 fprintf(stderr,"\nUnrecognized option %s\n",argv[i]);
451 MMG2D_usage(argv[0]);
452 return 0;
453 }
454
455 }
456
457 else {
458 if ( mesh->namein == NULL ) {
459 if ( !MMG2D_Set_inputMeshName(mesh,argv[i]) )
460 return 0;
461 if ( mesh->info.imprim == -99 ) {
463 return 0;
464 }
465 }
466 else if ( mesh->nameout == NULL ) {
467 if ( !MMG2D_Set_outputMeshName(mesh,argv[i]) )
468 return 0;
469 }
470 else {
471 fprintf(stdout,"\nArgument %s ignored\n",argv[i]);
472 MMG2D_usage(argv[0]);
473 return 0;
474 }
475 }
476 i++;
477 }
478
480 if ( mesh->info.imprim == -99 ) {
481 fprintf(stdout,"\n -- PRINT (0 10(advised) -10) ?\n");
482 fflush(stdin);
483 MMG_FSCANF(stdin,"%d",&i);
485 return 0;
486 }
487
488 if ( mesh->namein == NULL ) {
489 fprintf(stdout,"\n -- INPUT MESH NAME ?\n");
490 fflush(stdin);
491 MMG_FSCANF(stdin,"%127s",namein);
492 if ( !MMG2D_Set_inputMeshName(mesh,namein) )
493 return 0;
494 }
495 if ( mesh->nameout == NULL ) {
497 return 0;
498 }
499
500 /* adp mode: if the metric name has been stored in sol, move it in met */
501 if ( met->namein==NULL && sol && sol->namein && !(mesh->info.iso || mesh->info.isosurf || mesh->info.lag>=0) ) {
503 return 0;
505 }
506
507 /* default : store solution (resp. displacement) name in iso
508 * (resp. lagrangian) mode, metric name otherwise */
509 tmp = ( mesh->info.iso || mesh->info.isosurf || mesh->info.lag >=0 ) ? sol : met;
510 assert ( tmp );
511 if ( tmp->namein == NULL ) {
513 return 0;
514 }
515 if ( met->nameout == NULL ) {
516 if ( !MMG2D_Set_outputSolName(mesh,met,"") )
517 return 0;
518 }
519
520 return 1;
521}
522
531
533
534 fprintf(stdout,"\n\n");
535
536 return 1;
537}
538
549 int ret,i,j,npar,nbr,split;
550 MMG5_int ref,rin,rex,br;
551 float fp1,fp2,fp3;
552 char *ptr,data[256];
553 FILE *in;
554 fpos_t position;
555
556 /* Check for parameter file */
557 if (mesh->info.fparam) {
558 strcpy(data,mesh->info.fparam);
559 }
560 else {
561 strcpy(data,mesh->namein);
562 }
563
564 ptr = MMG5_Get_filenameExt(data);
565
566 if ( ptr ) *ptr = '\0';
567 strcat(data,".mmg2d");
568
569 in = fopen(data,"rb");
570
571 if ( !in ) {
572 if ( !mesh->info.fparam ) {
573 sprintf(data,"%s","DEFAULT.mmg2d");
574 in = fopen(data,"rb");
575 if ( !in )
576 return 1;
577 }
578 else if (mesh->info.fparam ) {
579 fprintf(stderr," ** In %s: %s file NOT FOUND. \n",__func__,data);
580 fprintf(stdout," ## ERROR: UNABLE TO LOAD PARAMETER FILE.\n");
581 return 0;
582 }
583 }
584 if ( mesh->info.imprim >= 0 ) {
585 fprintf(stdout,"\n %%%% %s OPENED\n",data);
586 }
587
588 /* Read parameters */
589 while ( !feof(in) ) {
590 ret = fscanf(in,"%255s",data);
591 if ( !ret || feof(in) ) break;
592 for (i=0; (size_t)i<strlen(data); i++) data[i] = tolower(data[i]);
593
594 /* Read user-defined references for the LS mode */
595 if ( !strcmp(data,"lsreferences") ) {
596 ret = fscanf(in,"%d",&npar);
597 if ( !ret ) {
598 fprintf(stderr," %%%% Wrong format for lsreferences: %d\n",npar);
599 return 0;
600 }
601
603 return 0;
604 }
605 for (i=0; i<mesh->info.nmat; i++) {
606 MMG_FSCANF(in,"%" MMG5_PRId "",&ref);
607 fgetpos(in,&position);
608 MMG_FSCANF(in,"%255s",data);
609 split = MMG5_MMAT_NoSplit;
610 rin = rex = ref;
611 if ( strcmp(data,"nosplit") ) {
612 fsetpos(in,&position);
613 split = MMG5_MMAT_Split;
614 MMG_FSCANF(in,"%" MMG5_PRId "",&rin);
615 MMG_FSCANF(in,"%" MMG5_PRId "",&rex);
616 }
617 if ( !MMG2D_Set_multiMat(mesh,met,ref,split,rin,rex) ) {
618 return 0;
619 }
620 }
621 }
622 /* Read user-defined local parameters and store them in the structure info->par */
623 else if ( !strcmp(data,"parameters") ) {
624 ret = fscanf(in,"%d",&npar);
625
626 if ( !ret ) {
627 fprintf(stderr," %%%% Wrong format for parameters: %d\n",npar);
628 return 0;
629 }
630 else if ( npar > MMG5_LPARMAX ) {
631 fprintf(stderr," %%%% Too many local parameters %d. Abort\n",npar);
632 return 0;
633 }
634
635 /* Allocate memory and fill the info->par table (adding one, corresponding to the command line data) */
636 if ( npar ) {
638 return 0;
639
640 for (i=0; i<mesh->info.npar; i++) {
641 ret = fscanf(in,"%" MMG5_PRId " %255s",&ref,data);
642 if ( ret ) ret = fscanf(in,"%f %f %f",&fp1,&fp2,&fp3);
643
644 if ( !ret ) {
645 fprintf(stderr," %%%% Wrong format: %s\n",data);
646 return (0);
647 }
648
649 for (j=0; (size_t)j<strlen(data); j++) data[j] = tolower(data[j]);
650 if ( !strcmp(data,"triangles") || !strcmp(data,"triangle") ) {
651 if ( !MMG2D_Set_localParameter(mesh,met,MMG5_Triangle,ref,fp1,fp2,fp3) ) {
652 return 0;
653 }
654 }
655 else if ( !strcmp(data,"edges") || !strcmp(data,"edge") ) {
656 if ( !MMG2D_Set_localParameter(mesh,met,MMG5_Edg,ref,fp1,fp2,fp3) ) {
657 return 0;
658 }
659 }
660 else {
661 fprintf(stderr," %%%% Wrong format: %s\n",data);
662 return 0;
663 }
664 }
665 }
666 }
667 /* Read user-defined references where connected components should stay attached in ls mode */
668 else if ( !strcmp(data,"lsbasereferences") ) {
669 MMG_FSCANF(in,"%d",&nbr);
671 return 0;
672
673 for (i=0; i<mesh->info.nbr; i++) {
674 MMG_FSCANF(in,"%" MMG5_PRId "",&br);
675 if ( !MMG2D_Set_lsBaseReference(mesh,met,br) ) {
676 return 0;
677 }
678 }
679 }
680 else {
681 fprintf(stderr," %%%% Wrong format: %s\n",data);
682 return 0;
683 }
684 }
685
686 fclose(in);
687 return 1;
688}
689
690/* Free the structure dedicated to the management of multiple local parameters */
692
693 free(mesh->info.par);
694 mesh->info.npar = 0;
695
696 return 1;
697}
698
700 MMG5_pTria pt,pt1;
701 MMG5_pEdge ped;
702 MMG5_int *adja,k,i,j,i1,i2,iel;
703
704 *nb_edges = 0;
705 if ( mesh->tria ) {
706 /* Create the triangle adjacency if needed */
707 if ( !mesh->adja ) {
708 if ( !MMG2D_hashTria( mesh ) ) {
709 fprintf(stderr,"\n ## Error: %s: unable to create "
710 "adjacency table.\n",__func__);
711 return 0;
712 }
713 }
714
715 /* Count the number of non boundary edges */
716 for ( k=1; k<=mesh->nt; k++ ) {
717 pt = &mesh->tria[k];
718 if ( !MG_EOK(pt) ) continue;
719
720 adja = &mesh->adja[3*(k-1)+1];
721
722 for ( i=0; i<3; i++ ) {
723 iel = adja[i] / 3;
724 assert ( iel != k );
725
726 pt1 = &mesh->tria[iel];
727
728 if ( (!iel) || (pt->ref != pt1->ref) ||
729 ((pt->ref==pt1->ref) && MG_SIN(pt->tag[i])) ||
730 (mesh->info.opnbdy && pt->tag[i]) ) {
731 /* Do not treat boundary edges */
732 continue;
733 }
734 if ( k < iel ) {
735 /* Treat edge from the triangle with lowest index */
736 ++(*nb_edges);
737 }
738 }
739 }
740
741 /* Append the non boundary edges to the boundary edges array */
742 if ( mesh->namax ) {
743 MMG5_ADD_MEM(mesh,(*nb_edges)*sizeof(MMG5_Edge),"non boundary edges",
744 printf(" Exit program.\n");
745 return 0);
746 MMG5_SAFE_RECALLOC(mesh->edge,(mesh->namax+1),(mesh->namax+(*nb_edges)+1),
747 MMG5_Edge,"non bdy edges arrray",return 0);
748 }
749 else {
750 MMG5_ADD_MEM(mesh,((*nb_edges)+1)*sizeof(MMG5_Edge),"non boundary edges",
751 printf(" Exit program.\n");
752 return 0);
753 MMG5_SAFE_RECALLOC(mesh->edge,0,((*nb_edges)+1),
754 MMG5_Edge,"non bdy edges arrray",return 0);
755 }
756
757 j = mesh->namax+1;
758 for ( k=1; k<=mesh->nt; k++ ) {
759 pt = &mesh->tria[k];
760 if ( !MG_EOK(pt) ) continue;
761
762 adja = &mesh->adja[3*(k-1)+1];
763
764 for ( i=0; i<3; i++ ) {
765 i1 = MMG5_inxt2[i];
766 i2 = MMG5_iprv2[i];
767 iel = adja[i] / 3;
768 assert ( iel != k );
769
770 pt1 = &mesh->tria[iel];
771
772 if ( (!iel) || (pt->ref != pt1->ref) ||
773 ((pt->ref==pt1->ref) && MG_SIN(pt->tag[i])) ||
774 (mesh->info.opnbdy && pt->tag[i]) ) {
775 /* Do not treat boundary edges */
776 continue;
777 }
778 if ( k < iel ) {
779 /* Treat edge from the triangle with lowest index */
780 ped = &mesh->edge[j++];
781 assert ( ped );
782 ped->a = pt->v[i1];
783 ped->b = pt->v[i2];
784 ped->ref = pt->edg[i];
785 }
786 }
787 }
788 }
789 return 1;
790}
791
792int MMG2D_Get_nonBdyEdge(MMG5_pMesh mesh, MMG5_int* e0, MMG5_int* e1, MMG5_int* ref, MMG5_int idx) {
793 MMG5_pEdge ped;
794 size_t na_tot=0;
795 char *ptr_c = (char*)mesh->edge;
796
797 if ( !mesh->edge ) {
798 fprintf(stderr,"\n ## Error: %s: edge array is not allocated.\n"
799 " Please, call the MMG2D_Get_numberOfNonBdyEdges function"
800 " before the %s one.\n",
801 __func__,__func__);
802 return 0;
803 }
804
805 ptr_c = ptr_c-sizeof(size_t);
806 na_tot = *((size_t*)ptr_c);
807
808 if ( mesh->namax==(MMG5_int)na_tot ) {
809 fprintf(stderr,"\n ## Error: %s: no internal edge.\n"
810 " Please, call the MMG2D_Get_numberOfNonBdyEdges function"
811 " before the %s one and check that the number of internal"
812 " edges is non null.\n",
813 __func__,__func__);
814 return 0;
815 }
816
817 if ( mesh->namax+idx > (MMG5_int)na_tot ) {
818 fprintf(stderr,"\n ## Error: %s: Can't get the internal edge of index %" MMG5_PRId "."
819 " Index must be between 1 and %"MMG5_PRId".\n",
820 __func__,idx,(MMG5_int)na_tot-mesh->namax);
821 return 0;
822 }
823
824 ped = &mesh->edge[mesh->namax+idx];
825
826 *e0 = ped->a;
827 *e1 = ped->b;
828
829 if ( ref != NULL ) {
830 *ref = mesh->edge[mesh->namax+idx].ref;
831 }
832
833 return 1;
834}
835
836int MMG2D_Get_adjaTri(MMG5_pMesh mesh, MMG5_int kel, MMG5_int listri[3]) {
837
838 if ( ! mesh->adja ) {
839 if (! MMG2D_hashTria(mesh))
840 return 0;
841 }
842
843 listri[0] = mesh->adja[3*(kel-1)+1]/3;
844 listri[1] = mesh->adja[3*(kel-1)+2]/3;
845 listri[2] = mesh->adja[3*(kel-1)+3]/3;
846
847 return 1;
848}
849
850MMG5_int MMG2D_Get_adjaVertices(MMG5_pMesh mesh, MMG5_int ip, MMG5_int lispoi[MMG2D_LMAX])
851{
852 MMG5_int start;
853
854 if ( !mesh->tria ) return 0;
855
856 start=MMG2D_findTria(mesh,ip);
857 if ( !start ) return 0;
858
859 return MMG2D_Get_adjaVerticesFast(mesh,ip,start,lispoi);
860}
861
862MMG5_int MMG2D_Get_adjaVerticesFast(MMG5_pMesh mesh, MMG5_int ip,MMG5_int start, MMG5_int lispoi[MMG2D_LMAX])
863{
864 MMG5_pTria pt;
865 int iploc,i,i1,i2;
866 MMG5_int prevk,k,*adja,nbpoi;
867
868 pt = &mesh->tria[start];
869
870 for ( iploc=0; iploc<3; ++iploc ) {
871 if ( pt->v[iploc] == ip ) break;
872 }
873
874 assert(iploc!=3);
875
876 k = start;
877 i = iploc;
878 nbpoi = 0;
879 do {
880 if ( nbpoi == MMG2D_LMAX ) {
881 fprintf(stderr,"\n ## Warning: %s: unable to compute adjacent"
882 " vertices of the vertex %" MMG5_PRId ":\nthe ball of point contain too many"
883 " elements.\n",__func__,ip);
884 return 0;
885 }
886 i1 = MMG5_inxt2[i];
887 lispoi[nbpoi] = mesh->tria[k].v[i1];
888 ++nbpoi;
889
890 adja = &mesh->adja[3*(k-1)+1];
891 prevk = k;
892 k = adja[i1] / 3;
893 i = adja[i1] % 3;
894 i = MMG5_inxt2[i];
895 }
896 while ( k && k != start );
897
898 if ( k > 0 ) return nbpoi;
899
900 /* store the last point of the boundary triangle */
901 if ( nbpoi == MMG2D_LMAX ) {
902 fprintf(stderr,"\n ## Warning: %s: unable to compute adjacent vertices of the"
903 " vertex %" MMG5_PRId ":\nthe ball of point contain too many elements.\n",
904 __func__,ip);
905 return 0;
906 }
907 i1 = MMG5_inxt2[i1];
908 lispoi[nbpoi] = mesh->tria[prevk].v[i1];
909 ++nbpoi;
910
911 /* check if boundary hit */
912 k = start;
913 i = iploc;
914 do {
915 adja = &mesh->adja[3*(k-1)+1];
916 i2 = MMG5_iprv2[i];
917 k = adja[i2] / 3;
918 if ( k == 0 ) break;
919
920 if ( nbpoi == MMG2D_LMAX ) {
921 fprintf(stderr,"\n ## Warning: %s: unable to compute adjacent vertices of the"
922 " vertex %" MMG5_PRId ":\nthe ball of point contain too many elements.\n",
923 __func__,ip);
924 return 0;
925 }
926 i = adja[i2] % 3;
927 lispoi[nbpoi] = mesh->tria[k].v[i];
928 ++nbpoi;
929
930 i = MMG5_iprv2[i];
931 }
932 while ( k );
933
934 return nbpoi;
935}
936
937int MMG2D_Get_triFromEdge(MMG5_pMesh mesh, MMG5_int ked, MMG5_int *ktri, int *ied)
938{
939 MMG5_int val;
940
941 val = mesh->edge[ked].base;
942
943 if ( !val ) {
944 fprintf(stderr," ## Error: %s: the main fonction of the Mmg library must be"
945 " called before this function.\n",__func__);
946 return 0;
947 }
948
949 *ktri = val/3;
950
951 *ied = val%3;
952
953 return 1;
954}
955
956int MMG2D_Get_trisFromEdge(MMG5_pMesh mesh, MMG5_int ked, MMG5_int ktri[2], int ied[2])
957{
958 int ier;
959 MMG5_int itri;
960#ifndef NDEBUG
961 MMG5_int ia0,ib0,ia1,ib1;
962#endif
963
964 ktri[0] = ktri[1] = 0;
965 ied[0] = ied[1] = 0;
966
967 ier = MMG2D_Get_triFromEdge(mesh, ked, ktri, ied);
968
969 if ( !ier ) return 0;
970
971 if ( !mesh->adja ) {
972 if (!MMG2D_hashTria(mesh) )
973 return 0;
974 }
975
976 itri = mesh->adja[3*(*ktri-1) + *ied + 1 ];
977
978 if ( itri ) {
979 ktri[1] = itri/3;
980 ied[1] = itri%3;
981
982#ifndef NDEBUG
983 ia0 = mesh->tria[ktri[0]].v[MMG5_inxt2[ied[0]]];
984 ib0 = mesh->tria[ktri[0]].v[MMG5_iprv2[ied[0]]];
985
986 ia1 = mesh->tria[ktri[1]].v[MMG5_inxt2[ied[1]]];
987 ib1 = mesh->tria[ktri[1]].v[MMG5_iprv2[ied[1]]];
988
989 assert ( ( (ia0 == ia1) && (ib0 == ib1) ) ||
990 ( (ia0 == ib1) && (ib0 == ia1) ) );
991#endif
992 }
993
994 return 1;
995}
996
998 double hsiz;
999
1000 /* Set solution size */
1001 if ( mesh->info.ani ) {
1002 met->size = 3;
1003 }
1004 else {
1005 met->size = 1;
1006 }
1007
1008 /* Memory alloc */
1009 if ( !MMG2D_Set_solSize(mesh,met,MMG5_Vertex,mesh->np,met->size) )
1010 return 0;
1011
1012 if ( !MMG5_Compute_constantSize(mesh,met,&hsiz) )
1013 return 0;
1014
1015 mesh->info.hsiz = hsiz;
1016
1017 MMG5_Set_constantSize(mesh,met,hsiz);
1018
1019 return 1;
1020}
1021
1022int MMG2D_Compute_eigenv(double m[3],double lambda[2],double vp[2][2]) {
1023
1024 return MMG5_eigensym(m,lambda,vp);
1025
1026}
1027
1028
1030 MMG5_int k;
1031
1032 for ( k=1; k<=mesh->np; ++k ) {
1033 mesh->point[k].tag = 0;
1034 }
1035
1036}
1037
1039
1040 if ( mesh->adja )
1042
1043 if ( mesh->tria )
1045
1046 mesh->nt = 0;
1047 mesh->nti = 0;
1048 mesh->nenil = 0;
1049
1050 return;
1051}
1052
1054
1055 if ( mesh->edge )
1057
1058 if ( mesh->xpoint )
1060
1061 mesh->na = 0;
1062 mesh->nai = 0;
1063 mesh->nanil = 0;
1064
1065 mesh->xp = 0;
1066
1067 return;
1068}
1069
1071
1072 /* sol */
1073 if ( !sol ) return;
1074
1075 if ( sol->m )
1076 MMG5_DEL_MEM(mesh,sol->m);
1077
1078 if ( sol->namein ) {
1080 }
1081
1082 if ( sol->nameout ) {
1084 }
1085
1086 memset ( sol, 0x0, sizeof(MMG5_Sol) );
1087
1088 /* Reset state to a scalar status */
1089 sol->dim = 2;
1090 sol->ver = 2;
1091 sol->size = 1;
1092 sol->type = 1;
1093
1094 return;
1095}
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 MMG2D_Set_lsBaseReference(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int br)
Set level-set base reference.
int MMG2D_Set_iparameter(MMG5_pMesh mesh, MMG5_pSol sol, int iparam, MMG5_int val)
Set integer parameter iparam to value val.
int MMG2D_Set_solSize(MMG5_pMesh mesh, MMG5_pSol sol, int typEntity, MMG5_int np, int typSol)
Set the size and type of a solution field.
int MMG2D_Set_inputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solin)
Set the name of the input solution file.
int MMG2D_Set_inputMeshName(MMG5_pMesh mesh, const char *meshin)
Set the name of the input mesh.
int MMG2D_Set_localParameter(MMG5_pMesh mesh, MMG5_pSol sol, int typ, MMG5_int ref, double hmin, double hmax, double hausd)
Set local parameters.
int MMG2D_Set_dparameter(MMG5_pMesh mesh, MMG5_pSol sol, int dparam, double val)
Set double parameter dparam to value val.
int MMG2D_Set_multiMat(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int ref, int split, MMG5_int rin, MMG5_int rout)
Set the reference mapping for the elements of ref ref in LS discretization mode.
int MMG2D_Set_inputParamName(MMG5_pMesh mesh, const char *fparamin)
Set the name of the input parameter file.
int MMG2D_Set_outputMeshName(MMG5_pMesh mesh, const char *meshout)
Set the name of the output mesh file.
int MMG2D_Set_outputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solout)
Set the name of the output solution file.
int ier
tmp[*strlen0]
MMG5_pMesh MMG5_pSol * sol
if(!ier) exit(EXIT_FAILURE)
MMG5_pMesh * mesh
int MMG5_gradsizreq_ani(MMG5_pMesh mesh, MMG5_pSol met)
Definition: anisosiz.c:2322
int MMG5_compute_meanMetricAtMarkedPoints_ani(MMG5_pMesh mesh, MMG5_pSol met)
Definition: anisosiz.c:2205
int MMG2D_defsiz_ani(MMG5_pMesh mesh, MMG5_pSol met)
Definition: anisosiz_2d.c:363
int MMG5_eigensym(double m[3], double lambda[2], double vp[2][2])
Definition: eigenv.c:973
int MMG2D_hashTria(MMG5_pMesh mesh)
Definition: hash_2d.c:35
int MMG2D_intmet_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int k, int8_t i, MMG5_int ip, double s)
Definition: intmet_2d.c:38
int MMG2D_intmet_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int k, int8_t i, MMG5_int ip, double s)
Definition: intmet_2d.c:209
int MMG5_gradsiz_iso(MMG5_pMesh mesh, MMG5_pSol met)
Definition: isosiz.c:278
int MMG5_compute_meanMetricAtMarkedPoints_iso(MMG5_pMesh mesh, MMG5_pSol met)
Definition: isosiz.c:167
int MMG5_gradsizreq_iso(MMG5_pMesh mesh, MMG5_pSol met)
Definition: isosiz.c:372
int MMG2D_defsiz_iso(MMG5_pMesh mesh, MMG5_pSol met)
Definition: isosiz_2d.c:133
double MMG2D_lencurv_iso(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int ip1, MMG5_int ip2)
Definition: length_2d.c:63
double MMG2D_lencurv_ani(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int ip1, MMG5_int ip2)
Definition: length_2d.c:82
API headers and documentation for the mmg2d library.
@ MMG2D_DPARAM_hgradreq
Definition: libmmg2d.h:140
@ MMG2D_IPARAM_numsubdomain
Definition: libmmg2d.h:128
@ MMG2D_IPARAM_iso
Definition: libmmg2d.h:116
@ MMG2D_IPARAM_optim
Definition: libmmg2d.h:121
@ MMG2D_IPARAM_numberOfLocalParam
Definition: libmmg2d.h:129
@ MMG2D_IPARAM_nosurf
Definition: libmmg2d.h:125
@ MMG2D_DPARAM_hgrad
Definition: libmmg2d.h:139
@ MMG2D_DPARAM_hmin
Definition: libmmg2d.h:135
@ MMG2D_IPARAM_mem
Definition: libmmg2d.h:113
@ MMG2D_DPARAM_xreg
Definition: libmmg2d.h:142
@ MMG2D_IPARAM_nofem
Definition: libmmg2d.h:144
@ MMG2D_IPARAM_nosizreq
Definition: libmmg2d.h:133
@ MMG2D_IPARAM_isoref
Definition: libmmg2d.h:145
@ MMG2D_DPARAM_rmc
Definition: libmmg2d.h:143
@ MMG2D_DPARAM_hausd
Definition: libmmg2d.h:138
@ MMG2D_IPARAM_angle
Definition: libmmg2d.h:115
@ MMG2D_DPARAM_hmax
Definition: libmmg2d.h:136
@ MMG2D_IPARAM_isosurf
Definition: libmmg2d.h:117
@ MMG2D_DPARAM_ls
Definition: libmmg2d.h:141
@ MMG2D_IPARAM_noinsert
Definition: libmmg2d.h:122
@ MMG2D_IPARAM_nreg
Definition: libmmg2d.h:126
@ MMG2D_IPARAM_noswap
Definition: libmmg2d.h:123
@ MMG2D_IPARAM_xreg
Definition: libmmg2d.h:127
@ MMG2D_IPARAM_nomove
Definition: libmmg2d.h:124
@ MMG2D_IPARAM_lag
Definition: libmmg2d.h:119
@ MMG2D_IPARAM_opnbdy
Definition: libmmg2d.h:118
@ MMG2D_IPARAM_verbose
Definition: libmmg2d.h:112
@ MMG2D_IPARAM_3dMedit
Definition: libmmg2d.h:120
@ MMG2D_IPARAM_numberOfLSBaseReferences
Definition: libmmg2d.h:130
@ MMG2D_IPARAM_numberOfMat
Definition: libmmg2d.h:131
@ MMG2D_DPARAM_angleDetection
Definition: libmmg2d.h:134
@ MMG2D_DPARAM_hsiz
Definition: libmmg2d.h:137
@ MMG2D_IPARAM_debug
Definition: libmmg2d.h:114
LIBMMG2D_EXPORT int(* MMG2D_doSol)(MMG5_pMesh mesh, MMG5_pSol met)
Compute unit tensor according to the lengths of the edges passing through a vertex.
Definition: mmg2dexterns.c:9
#define MMG2D_LMAX
Definition: libmmg2d.h:98
MMG5_int MMG2D_findTria(MMG5_pMesh, MMG5_int)
Definition: locate_2d.c:221
int lissmet_ani(MMG5_pMesh mesh, MMG5_pSol sol)
Definition: lissmet_2d.c:47
double MMG2D_caltri_ani(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_pTria)
Definition: quality_2d.c:121
int MMG2D_doSol_iso(MMG5_pMesh mesh, MMG5_pSol sol)
Definition: solmap_2d.c:97
int MMG2D_doSol_ani(MMG5_pMesh mesh, MMG5_pSol sol)
Definition: solmap_2d.c:175
double MMG2D_caltri_iso(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_pTria)
Definition: quality_2d.c:107
int MMG2D_defaultValues(MMG5_pMesh mesh)
Print the default parameters values.
int MMG2D_usage(char *name)
Print help for mmg2d options.
int MMG2D_Get_trisFromEdge(MMG5_pMesh mesh, MMG5_int ked, MMG5_int ktri[2], int ied[2])
Find two triangles given the edge that they share.
int MMG2D_freeLocalPar(MMG5_pMesh mesh)
int MMG2D_Get_numberOfNonBdyEdges(MMG5_pMesh mesh, MMG5_int *nb_edges)
Get the number of non-boundary edges.
int MMG2D_Compute_eigenv(double m[3], double lambda[2], double vp[2][2])
Compute the real eigenvalues and eigenvectors of a symmetric matrix.
void MMG2D_setfunc(MMG5_pMesh mesh, MMG5_pSol met)
Set function pointers for length, caltri... depending if case is iso or aniso.
void MMG2D_Reset_verticestags(MMG5_pMesh mesh)
Reset the vertex tags.
MMG5_int MMG2D_Get_adjaVertices(MMG5_pMesh mesh, MMG5_int ip, MMG5_int lispoi[MMG2D_LMAX])
Return adjacent vertices of a triangle.
int MMG2D_Get_adjaTri(MMG5_pMesh mesh, MMG5_int kel, MMG5_int listri[3])
Return adjacent elements of a triangle.
MMG5_int MMG2D_Get_adjaVerticesFast(MMG5_pMesh mesh, MMG5_int ip, MMG5_int start, MMG5_int lispoi[MMG2D_LMAX])
Return adjacent vertices of a triangle.
void MMG2D_Free_edges(MMG5_pMesh mesh)
Free the mesh edges (and the associated xpoints).
int MMG2D_Get_nonBdyEdge(MMG5_pMesh mesh, MMG5_int *e0, MMG5_int *e1, MMG5_int *ref, MMG5_int idx)
Get vertices and reference of a non-boundary edge.
int MMG2D_parsop(MMG5_pMesh mesh, MMG5_pSol met)
Read a file containing Local parameters (.mmg2d extension)
int MMG2D_parsar(int argc, char *argv[], MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol)
Store command line arguments.
void MMG2D_Free_triangles(MMG5_pMesh mesh)
Free the mesh elements (and the adjacency information).
void MMG2D_Free_solutions(MMG5_pMesh mesh, MMG5_pSol sol)
Free the solution.
int MMG2D_Get_triFromEdge(MMG5_pMesh mesh, MMG5_int ked, MMG5_int *ktri, int *ied)
Find a triangle given an adjacent triangle and an edge number.
int MMG2D_Set_constantSize(MMG5_pMesh mesh, MMG5_pSol met)
Compute a constant size map according to the hsiz, hmin and hmax parameters.
#define MMG5_MMAT_Split
Definition: libmmgtypes.h:211
@ MMG5_Tensor
Definition: libmmgtypes.h:221
#define MMG5_MMAT_NoSplit
Definition: libmmgtypes.h:203
@ MMG5_Vertex
Definition: libmmgtypes.h:230
@ MMG5_Edg
Definition: libmmgtypes.h:231
@ 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 MG_EOK(pt)
#define MMG5_ADD_MEM(mesh, size, message, law)
static const uint8_t MMG5_iprv2[3]
#define MG_SIN(tag)
#define MMG5_UNSET
static const uint8_t MMG5_inxt2[6]
#define MMG5_LPARMAX
#define MMG5_FILESTR_LGTH
#define MMG_FSCANF(stream, format,...)
#define MMG5_DEL_MEM(mesh, ptr)
#define MMG5_SAFE_RECALLOC(ptr, prevSize, newSize, type, message, law)
Structure to store edges of am MMG mesh.
Definition: libmmgtypes.h:311
MMG5_int b
Definition: libmmgtypes.h:312
MMG5_int base
Definition: libmmgtypes.h:314
MMG5_int ref
Definition: libmmgtypes.h:313
MMG5_int a
Definition: libmmgtypes.h:312
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
int8_t lag
Definition: libmmgtypes.h:547
char * fparam
Definition: libmmgtypes.h:555
MMG5_pPar par
Definition: libmmgtypes.h:524
MMG mesh structure.
Definition: libmmgtypes.h:613
MMG5_Info info
Definition: libmmgtypes.h:659
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_pxPoint xpoint
Definition: libmmgtypes.h:650
MMG5_int nenil
Definition: libmmgtypes.h:630
MMG5_int xp
Definition: libmmgtypes.h:628
MMG5_int namax
Definition: libmmgtypes.h:620
MMG5_int nt
Definition: libmmgtypes.h:620
MMG5_pTria tria
Definition: libmmgtypes.h:655
MMG5_int np
Definition: libmmgtypes.h:620
MMG5_pEdge edge
Definition: libmmgtypes.h:657
MMG5_int nanil
Definition: libmmgtypes.h:631
MMG5_int nti
Definition: libmmgtypes.h:620
MMG5_int nai
Definition: libmmgtypes.h:620
MMG5_int na
Definition: libmmgtypes.h:620
char * namein
Definition: libmmgtypes.h:660
uint16_t tag
Definition: libmmgtypes.h:290
char * nameout
Definition: libmmgtypes.h:683
char * namein
Definition: libmmgtypes.h:682
Structure to store triangles of a MMG mesh.
Definition: libmmgtypes.h:338
MMG5_int edg[3]
Definition: libmmgtypes.h:345
MMG5_int ref
Definition: libmmgtypes.h:341
uint16_t tag[3]
Definition: libmmgtypes.h:348
MMG5_int v[3]
Definition: libmmgtypes.h:340