Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
API_functions_2d.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
40#include "libmmg2d.h"
41#include "libmmg2d_private.h"
42
43int MMG2D_Init_mesh(const int starter,...) {
44 va_list argptr;
45 int ier;
46
48
50
52
53 return ier;
54}
55
62
67
69 return MMG5_Set_inputSolName(mesh,sol,solin);
70}
71
72int MMG2D_Set_outputMeshName(MMG5_pMesh mesh, const char* meshout) {
73
74 return MMG5_Set_outputMeshName(mesh,meshout);
75}
76
78 return MMG5_Set_outputSolName(mesh,sol,solout);
79}
80
82
83 /* Init common parameters for mmg2d, mmgs and mmg2d. */
85
86 /* default values for integers */
90 /* [0/1] ,avoid/allow surface modifications */
92 /* [0/1/2] , set 3D mode for 2D mesh */
94 /* [0/1] , set off/on normal regularization */
96 /* [0/1] , set off/on vertex regularization */
98 /* default values for doubles */
99 /* level set value */
100 mesh->info.ls = MMG5_LS;
101
102 /* Ridge detection */
104}
105
106
107int MMG2D_Set_iparameter(MMG5_pMesh mesh, MMG5_pSol sol, int iparam, MMG5_int val){
108 int k;
109
110 switch ( iparam ) {
111 /* Integer parameters */
113 mesh->info.imprim = val;
114 break;
115 case MMG2D_IPARAM_mem :
116 if ( val <= 0 ) {
117 fprintf(stderr,"\n ## Warning: %s: maximal memory authorized must"
118 " be strictly positive.\n",__func__);
119 fprintf(stderr," Reset to default value.\n");
120 }
121 else
122 mesh->info.mem = val;
123 if ( !MMG2D_memOption(mesh) ) return 0;
124 break;
125 case MMG2D_IPARAM_debug :
126 mesh->info.ddebug = val;
127 break;
128 case MMG2D_IPARAM_angle :
129 /* free table that may contains old ridges */
130 if ( mesh->htab.geom )
132 if ( mesh->xpoint )
134 if ( mesh->xtetra )
136 if ( !val )
137 mesh->info.dhd = MMG5_NR;
138 else {
139 if ( (mesh->info.imprim > 5) || mesh->info.ddebug )
140 fprintf(stderr,"\n ## Warning: %s: angle detection parameter"
141 " set to default value\n",__func__);
143 }
144 break;
145 case MMG2D_IPARAM_nofem :
146 mesh->info.setfem = (val==1)? 0 : 1;
147 break;
149 mesh->info.opnbdy = val;
150 break;
151 case MMG2D_IPARAM_iso :
152 mesh->info.iso = val;
153 break;
155 mesh->info.isoref = val;
156 break;
158 mesh->info.isosurf = val;
159 break;
160 case MMG2D_IPARAM_lag :
161#ifdef USE_ELAS
162 if ( val < 0 || val > 2 )
163 return 0;
164 mesh->info.lag = val;
165 /* No connectivity changes unless lag >= 2 */
166 if ( val < 2 ) {
168 return 0;
169 }
170#else
171 fprintf(stderr,"\n ## Error: %s"
172 " \"lagrangian motion\" option unavailable (-lag):\n"
173 " set the USE_ELAS CMake's flag to ON when compiling the mmg3d"
174 " library to enable this feature.\n",__func__);
175 return 0;
176#endif
177 break;
179 mesh->info.renum = val;
180 break;
182 mesh->info.nsd = val;
183 break;
184 case MMG2D_IPARAM_optim :
185 mesh->info.optim = val;
186 break;
188 mesh->info.noinsert = val;
189 break;
191 mesh->info.noswap = val;
192 break;
194 mesh->info.nomove = val;
195 break;
197 mesh->info.nosurf = val;
198 break;
199 case MMG2D_IPARAM_nreg :
200 mesh->info.nreg = val;
201 break;
202 case MMG2D_IPARAM_xreg :
203 mesh->info.xreg = val;
204 break;
206 mesh->info.nosizreq = val;
207 break;
209 if ( mesh->info.par ) {
211 if ( (mesh->info.imprim > 5) || mesh->info.ddebug )
212 fprintf(stderr,"\n ## Warning: %s: new local parameter values\n",__func__);
213 }
214 mesh->info.npar = val;
215 mesh->info.npari = 0;
216 mesh->info.parTyp = 0;
217
218 MMG5_ADD_MEM(mesh,mesh->info.npar*sizeof(MMG5_Par),"parameters",
219 printf(" Exit program.\n");
220 return 0);
222
223 MMG5_int inival;
224 if ( sizeof(MMG5_int) == 8 ) {
225 inival = LONG_MAX;
226 }
227 else {
228 inival = INT_MAX;
229 }
230
231 for (k=0; k<mesh->info.npar; k++) {
233 mesh->info.par[k].ref = inival;
235 mesh->info.par[k].hmin = mesh->info.hmin;
236 mesh->info.par[k].hmax = mesh->info.hmax;
237 }
238 break;
240 if ( mesh->info.br ) {
242 if ( (mesh->info.imprim > 5) || mesh->info.ddebug )
243 fprintf(stderr,"\n ## Warning: %s: new level-set based references values\n",__func__);
244 }
245 mesh->info.nbr = val;
246 mesh->info.nbri = 0;
247 MMG5_ADD_MEM(mesh,mesh->info.nbr*sizeof(MMG5_int),"References",
248 printf(" Exit program.\n");
249 return 0);
250 MMG5_SAFE_CALLOC(mesh->info.br,mesh->info.nbr,MMG5_int,return 0);
251
252 for (k=0; k<mesh->info.nbr; k++)
253 mesh->info.br[k] = 0;
254
255 break;
256
258 if ( mesh->info.mat ) {
260 if ( (mesh->info.imprim > 5) || mesh->info.ddebug )
261 fprintf(stderr,"\n ## Warning: %s: new multi materials values\n",__func__);
262 }
263 mesh->info.nmat = val;
264 mesh->info.nmati = 0;
265
266 MMG5_ADD_MEM(mesh,(mesh->info.nmat)*sizeof(MMG5_Mat),"multi material",
267 printf(" Exit program.\n");
268 return 0);
270 break;
271
273 mesh->info.ani = val;
274 break;
275 default :
276 fprintf(stderr,"\n ## Error: %s: unknown type of parameter\n",__func__);
277 return 0;
278 }
279 /* other options */
280
281 return 1;
282}
283
284int MMG2D_Set_dparameter(MMG5_pMesh mesh, MMG5_pSol sol, int dparam, double val){
285
286 switch ( dparam ) {
287 /* double parameters */
289 mesh->info.dhd = val;
290 mesh->info.dhd = MG_MAX(0.0, MG_MIN(180.0,mesh->info.dhd));
291 mesh->info.dhd = cos(mesh->info.dhd*M_PI/180.0);
292 break;
293 case MMG2D_DPARAM_hmin :
294 mesh->info.sethmin = 1;
295 mesh->info.hmin = val;
296 if ( mesh->info.sethmax && ( mesh->info.hmin >= mesh->info.hmax ) ) {
297 fprintf(stderr,"\n ## Error: hmin value must be strictly lower than hmax one"
298 " (hmin = %lf hmax = %lf ).\n",mesh->info.hmin, mesh->info.hmax);
299 return 0;
300 }
301 if ( val <= 0. ) {
302 fprintf(stderr,"\n ## Error: hmin must be strictly positive "
303 "(minimal edge length).\n");
304 return 0;
305 }
306
307 break;
308 case MMG2D_DPARAM_hmax :
309 mesh->info.sethmax = 1;
310 mesh->info.hmax = val;
311 if ( mesh->info.sethmin && ( mesh->info.hmin >= mesh->info.hmax ) ) {
312 fprintf(stderr,"\n ## Error: hmin value must be strictly lower than hmax one"
313 " (hmin = %lf hmax = %lf ).\n",mesh->info.hmin, mesh->info.hmax);
314 return 0;
315 }
316 if ( val <= 0. ) {
317 fprintf(stderr,"\n ## Error: hmax must be strictly positive "
318 "(maximal edge length).\n");
319 return 0;
320 }
321
322 break;
323 case MMG2D_DPARAM_hsiz :
324 mesh->info.hsiz = val;
325 break;
326 case MMG2D_DPARAM_hgrad :
327 mesh->info.hgrad = val;
328 if ( mesh->info.hgrad <= 0.0 ) {
330 }
331 else {
332 mesh->info.hgrad = log(mesh->info.hgrad);
333 }
334 break;
336 mesh->info.hgradreq = val;
337 if ( mesh->info.hgradreq <= 0.0 ) {
339 }
340 else {
342 }
343 break;
344 case MMG2D_DPARAM_hausd :
345 if ( val <=0 ) {
346 fprintf(stderr,"\n ## Error: %s: hausdorff number must be"
347 " strictly positive.\n",__func__);
348 return 0;
349 }
350 else
351 mesh->info.hausd = val;
352 break;
353 case MMG2D_DPARAM_ls :
354 mesh->info.ls = val;
355 break;
356 case MMG2D_DPARAM_rmc :
357 if ( !val ) {
358 /* Default value */
360 }
361 else {
362 /* User customized value */
363 mesh->info.rmc = val;
364 }
365 break;
366 default :
367 fprintf(stderr,"\n ## Error: %s: unknown type of parameter\n",
368 __func__);
369 return 0;
370 }
371 return 1;
372}
373
375 double hmin,double hmax,double hausd){
376 MMG5_pPar par;
377 int k;
378
379 if ( !mesh->info.npar ) {
380 fprintf(stderr,"\n ## Error: %s: You must set the number of local"
381 " parameters",__func__);
382 fprintf(stderr," with the MMG2D_Set_iparameters function before setting");
383 fprintf(stderr," values in local parameters structure. \n");
384 return 0;
385 }
386 if ( mesh->info.npari >= mesh->info.npar ) {
387 fprintf(stderr,"\n ## Error: %s: unable to set a new local parameter.\n",
388 __func__);
389 fprintf(stderr," max number of local parameters: %d\n",mesh->info.npar);
390 return 0;
391 }
392 if ( typ != MMG5_Triangle && typ != MMG5_Edg ) {
393 fprintf(stderr,"\n ## Warning: %s: you must apply your local parameters",
394 __func__);
395 fprintf(stderr," on triangles (MMG5_Triangle or %d) or edges"
396 " (MMG5_Edg or %d).\n",MMG5_Triangle,MMG5_Edg);
397 fprintf(stderr,"\n ## Unknown type of entity: ignored.\n");
398 return 0;
399 }
400 if ( ref < 0 ) {
401 fprintf(stderr,"\n ## Error: %s: negative references are not allowed.\n",
402 __func__);
403 return 0;
404 }
405 if ( hmin <= 0 ) {
406 fprintf(stderr,"\n ## Error: %s: negative hmin value is not allowed.\n",
407 __func__);
408 return 0;
409 }
410 if ( hmax <= 0 ) {
411 fprintf(stderr,"\n ## Error: %s: negative hmax value is not allowed.\n",
412 __func__);
413 return 0;
414 }
415 if ( hausd <= 0 ) {
416 fprintf(stderr,"\n ## Error: %s: negative hausd value is not allowed.\n",
417 __func__);
418 return 0;
419 }
420
421 for (k=0; k<mesh->info.npari; k++) {
422 par = &mesh->info.par[k];
423
424 if ( par->elt == typ && par->ref == ref ) {
425 par->hausd = hausd;
426 par->hmin = hmin;
427 par->hmax = hmax;
428 if ( (mesh->info.imprim > 5) || mesh->info.ddebug ) {
429 fprintf(stderr,"\n ## Warning: %s: new parameters (hausd, hmin and hmax)",
430 __func__);
431 fprintf(stderr," for entities of type %d and of ref %" MMG5_PRId "\n",typ,ref);
432 }
433 return 1;
434 }
435 }
436
437 mesh->info.par[mesh->info.npari].elt = typ;
438 mesh->info.par[mesh->info.npari].ref = ref;
439 mesh->info.par[mesh->info.npari].hmin = hmin;
440 mesh->info.par[mesh->info.npari].hmax = hmax;
441 mesh->info.par[mesh->info.npari].hausd = hausd;
442
443 switch ( typ )
444 {
445 case ( MMG5_Triangle ):
447 break;
448 case ( MMG5_Edg ):
450 break;
451 default:
452 fprintf(stderr,"\n ## Error: %s: unexpected entity type: %s.\n",
453 __func__,MMG5_Get_entitiesName(typ));
454 return 0;
455 }
456
457 mesh->info.npari++;
458
459 return 1;
460}
461
463 int split,MMG5_int rin,MMG5_int rout){
464 return MMG5_Set_multiMat(mesh,sol,ref,split,rin,rout);
465}
466
470
471
472int MMG2D_Set_meshSize(MMG5_pMesh mesh, MMG5_int np, MMG5_int nt, MMG5_int nquad, MMG5_int na) {
473
474 if ( ( (mesh->info.imprim > 5) || mesh->info.ddebug ) &&
475 ( mesh->point || mesh->tria || mesh->edge) )
476 fprintf(stderr,"\n ## Warning: %s: old mesh deletion.\n",__func__);
477
478 if ( mesh->point )
480 if ( mesh->tria )
482 if ( mesh->quadra )
484 if ( mesh->edge )
486
487 mesh->np = np;
488 mesh->nt = nt;
489 mesh->nquad = nquad;
490 mesh->na = na;
491 mesh->npi = mesh->np;
492 mesh->nti = mesh->nt;
493 mesh->nai = mesh->na;
494
495 /*tester si -m definie : renvoie 0 si pas ok et met la taille min dans info.mem */
496 if( mesh->info.mem > 0) {
497 if((mesh->npmax < mesh->np || mesh->ntmax < mesh->nt || mesh->namax < mesh->na) ) {
498 if ( !MMG2D_memOption(mesh) ) return 0;
499 } else if(mesh->info.mem < 39) {
500 fprintf(stderr,"\n ## Error: %s: not enough memory (%d).\n",
501 __func__,mesh->info.mem);
502 return 0;
503 }
504 } else {
505 if ( !MMG2D_memOption(mesh) ) return 0;
506 }
507
508 /* Mesh allocation and linkage */
509 if ( !MMG2D_setMeshSize_alloc( mesh ) ) return 0;
510
511 return 1;
512}
513
514int MMG2D_Set_solSize(MMG5_pMesh mesh, MMG5_pSol sol, int typEntity, MMG5_int np, int typSol) {
515
516 if ( ( (mesh->info.imprim > 5) || mesh->info.ddebug ) && sol->m )
517 fprintf(stderr,"\n ## Warning: %s: old solution deletion.\n",__func__);
518
519 if ( typEntity != MMG5_Vertex ) {
520 fprintf(stderr,"\n ## Error: %s: mmg2d need a solution imposed on vertices.\n",
521 __func__);
522 return 0;
523 }
524
525 sol->type = typSol;
526
527 if ( typSol == MMG5_Scalar ) {
528 sol->size = 1;
529 }
530 else if ( typSol == MMG5_Vector ) {
531 sol->size = 2;
532 }
533 else if ( typSol == MMG5_Tensor ) {
534 sol->size = 3;
535 }
536 else {
537 fprintf(stderr,"\n ## Error: %s: type of solution not yet implemented.\n",
538 __func__);
539 return 0;
540 }
541
542 sol->dim = 2;
543 if ( np ) {
544 sol->np = np;
545 sol->npi = np;
546 if ( sol->m )
548
549 sol->npmax = mesh->npmax;
550 MMG5_ADD_MEM(mesh,(sol->size*(sol->npmax+1))*sizeof(double),"initial solution",
551 printf(" Exit program.\n");
552 return 0);
553 MMG5_SAFE_CALLOC(sol->m,(sol->size*(sol->npmax+1)),double,return 0);
554 }
555 return 1;
556}
557
559 MMG5_int np, int *typSol) {
560 MMG5_pSol psl;
561 char data[18];
562 int j;
563
564 if ( ( (mesh->info.imprim > 5) || mesh->info.ddebug ) && mesh->nsols ) {
565 if ( *sol ) {
566 fprintf(stderr,"\n ## Warning: %s: old solutions array deletion.\n",
567 __func__);
569 }
570 }
571
573 mesh->nsols = nsols;
574
575 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
576 return 0);
577 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return 0);
578
579 for ( j=0; j<nsols; ++j ) {
580 psl = *sol + j;
581 psl->ver = 2;
582
583 /* Give an arbitrary name to the solution */
584 sprintf(data,"sol_%d",j);
585 if ( !MMG2D_Set_inputSolName(mesh,psl,data) ) {
586 return 0;
587 }
588 /* Give an arbitrary name to the solution */
589 sprintf(data,"sol_%d.o",j);
590 if ( !MMG2D_Set_outputSolName(mesh,psl,data) ) {
591 return 0;
592 }
593
594 if ( !MMG2D_Set_solSize(mesh,psl,MMG5_Vertex,mesh->np,typSol[j]) ) {
595 fprintf(stderr,"\n ## Error: %s: unable to set the size of the"
596 " solution num %d.\n",__func__,j);
597 return 0;
598 }
599 }
600 return 1;
601}
602
603int MMG2D_Get_solSize(MMG5_pMesh mesh, MMG5_pSol sol, int* typEntity, MMG5_int* np,
604 int* typSol) {
605
606 if ( typEntity != NULL )
607 *typEntity = MMG5_Vertex;
608
609 if ( typSol != NULL ) {
610 if ( sol->size == 1 )
611 *typSol = MMG5_Scalar;
612 else if ( sol->size == 2 )
613 *typSol = MMG5_Vector;
614 else if ( sol->size == 3 )
615 *typSol = MMG5_Tensor;
616 else
617 *typSol = MMG5_Notype;
618 }
619
620 assert( (!sol->np) || (sol->np == mesh->np));
621
622 if ( np != NULL )
623 *np = sol->np;
624
625 return 1;
626}
627
629 MMG5_int* np, int* typSol) {
630 MMG5_pSol psl;
631 MMG5_int j;
632
633 if ( !mesh ) {
634 fprintf(stderr,"\n ## Error: %s: your mesh structure must be allocated"
635 " and filled\n",__func__);
636 return 0;
637 }
638
639 if ( nsols != NULL )
640 *nsols = mesh->nsols;
641
642 for ( j=0; j<mesh->nsols; ++j ) {
643 psl = *sol + j;
644
645 if ( typSol != NULL ) {
646 typSol[j] = psl->type;
647 }
648
649 assert( (!psl->np) || (psl->np == mesh->np));
650 }
651 if ( np != NULL )
652 *np = mesh->np;
653
654 return 1;
655}
656
657int MMG2D_Get_meshSize(MMG5_pMesh mesh, MMG5_int* np, MMG5_int* nt, MMG5_int* nquad, MMG5_int* na) {
658 MMG5_int k;
659
660 if ( np != NULL )
661 *np = mesh->np;
662 if ( nt != NULL )
663 *nt = mesh->nt;
664
665 if ( nquad != NULL )
666 *nquad = mesh->nquad;
667
668 if ( na != NULL ) {
669 // Edges are not packed, thus we must count it.
670 *na = 0;
671 if ( mesh->na ) {
672 for (k=1; k<=mesh->na; k++) {
673 if ( mesh->edge[k].a ) ++(*na);
674 }
675 }
676 }
677
678 return 1;
679}
680
681int MMG2D_Set_vertex(MMG5_pMesh mesh, double c0, double c1, MMG5_int ref, MMG5_int pos) {
682
683 if ( !mesh->np ) {
684 fprintf(stderr,"\n ## Error: %s: you must set the number of points with the",
685 __func__);
686 fprintf(stderr," MMG2D_Set_meshSize function before setting vertices in mesh\n");
687 return 0;
688 }
689
690 if ( pos > mesh->npmax ) {
691 fprintf(stderr,"\n ## Error: %s: unable to allocate a new point.\n",
692 __func__);
693 fprintf(stderr," max number of points: %" MMG5_PRId "\n",mesh->npmax);
695 return 0;
696 }
697
698 if ( pos > mesh->np ) {
699 fprintf(stderr,"\n ## Error: %s: attempt to set new vertex at position %" MMG5_PRId ".",
700 __func__,pos);
701 fprintf(stderr," Overflow of the given number of vertices: %" MMG5_PRId "\n",mesh->np);
702 fprintf(stderr," ## Check the mesh size, its compactness or the position");
703 fprintf(stderr," of the vertex.\n");
704 return 0;
705 }
706
707 mesh->point[pos].c[0] = c0;
708 mesh->point[pos].c[1] = c1;
709 mesh->point[pos].ref = ref;
710 if ( mesh->nt )
711 mesh->point[pos].tag = MG_NUL;
712 else
713 mesh->point[pos].tag &= ~MG_NUL;
714
715 mesh->point[pos].flag = 0;
716 mesh->point[pos].tmp = 0;
717
718 return 1;
719}
720
722 assert ( k <= mesh->np );
723 mesh->point[k].tag |= MG_CRN;
724 return 1;
725}
726
728 assert ( k <= mesh->np );
729 mesh->point[k].tag &= ~MG_CRN;
730 return 1;
731}
732
734 assert ( k <= mesh->np );
735 mesh->point[k].tag |= MG_REQ;
736 mesh->point[k].tag &= ~MG_NUL;
737 return 1;
738}
739
741 assert ( k <= mesh->np );
742 mesh->point[k].tag &= ~MG_REQ;
743 return 1;
744}
745
746int MMG2D_Get_vertex(MMG5_pMesh mesh, double* c0, double* c1, MMG5_int* ref,
747 int* isCorner, int* isRequired) {
748
749 if ( mesh->npi == mesh->np ) {
750 mesh->npi = 0;
751 if ( mesh->info.ddebug ) {
752 fprintf(stderr,"\n ## Warning: %s: reset the internal counter of points.\n",
753 __func__);
754 fprintf(stderr," You must pass here exactly one time (the first time ");
755 fprintf(stderr,"you call the MMG2D_Get_vertex function).\n");
756 fprintf(stderr," If not, the number of call of this function");
757 fprintf(stderr," exceed the number of points: %" MMG5_PRId "\n ",mesh->np);
758 }
759 }
760
761 mesh->npi++;
762
763 if ( mesh->npi > mesh->np ) {
764 fprintf(stderr," ## Error: %s: unable to get point.\n",__func__);
765 fprintf(stderr," The number of call of MMG2D_Get_vertex function");
766 fprintf(stderr," exceed the number of points: %" MMG5_PRId "\n ",mesh->np);
767 return 0;
768 }
769
770 return MMG2D_GetByIdx_vertex( mesh,c0,c1,ref,isCorner,isRequired,mesh->npi);
771}
772
773int MMG2D_GetByIdx_vertex(MMG5_pMesh mesh, double* c0, double* c1, MMG5_int* ref,
774 int* isCorner, int* isRequired, MMG5_int idx) {
775
776 if ( idx < 1 || idx > mesh->np ) {
777 fprintf(stderr,"\n ## Error: %s: unable to get point at position %" MMG5_PRId ".\n",
778 __func__,idx);
779 fprintf(stderr," Your vertices numbering goes from 1 to %" MMG5_PRId "\n",mesh->np);
780 return 0;
781 }
782
783 *c0 = mesh->point[idx].c[0];
784 *c1 = mesh->point[idx].c[1];
785 if ( ref != NULL )
786 *ref = mesh->point[idx].ref;
787
788 if ( isCorner != NULL ) {
789 if ( mesh->point[idx].tag & MG_CRN )
790 *isCorner = 1;
791 else
792 *isCorner = 0;
793 }
794
795 if ( isRequired != NULL ) {
796 if ( mesh->point[idx].tag & MG_REQ )
797 *isRequired = 1;
798 else
799 *isRequired = 0;
800 }
801
802 return 1;
803}
804
805int MMG2D_Set_vertices(MMG5_pMesh mesh, double *vertices,MMG5_int *refs) {
806 MMG5_pPoint ppt;
807 MMG5_int i,j;
808
809 /*coordinates vertices*/
810 for (i=1;i<=mesh->np;i++)
811 {
812 ppt = &mesh->point[i];
813
814 j = (i-1)*2;
815 ppt->c[0] = vertices[j];
816 ppt->c[1] = vertices[j+1];
817
818 ppt->flag = 0;
819 ppt->tmp = 0;
820
821 if ( refs != NULL )
822 ppt->ref = refs[i-1];
823
824 if ( mesh->nt )
825 ppt->tag = MG_NUL;
826 else
827 ppt->tag &= ~MG_NUL;
828 }
829
830 return 1;
831}
832
833int MMG2D_Get_vertices(MMG5_pMesh mesh, double* vertices, MMG5_int* refs,
834 int* areCorners, int* areRequired) {
835 MMG5_pPoint ppt;
836 MMG5_int i,j;
837
838 for (i=1;i<=mesh->np;i++)
839 {
840 ppt = &mesh->point[i];
841
842 j = (i-1)*2;
843 vertices[j] = ppt->c[0];
844 vertices[j+1] = ppt->c[1];
845
846 j = i-1;
847 if ( refs != NULL )
848 refs[j] = ppt->ref;
849
850 if ( areCorners !=NULL ) {
851 if ( ppt->tag & MG_CRN )
852 areCorners[j] = 1;
853 else
854 areCorners[j] = 0;
855 }
856
857 if ( areRequired != NULL ) {
858 if ( ppt->tag & MG_REQ )
859 areRequired[j] = 1;
860 else
861 areRequired[j] = 0;
862 }
863 }
864
865 return 1;
866}
867
868int MMG2D_Set_triangle(MMG5_pMesh mesh, MMG5_int v0, MMG5_int v1, MMG5_int v2, MMG5_int ref, MMG5_int pos) {
869 MMG5_pPoint ppt;
870 MMG5_pTria pt;
871 double vol;
872 int i,j,ip;
873 MMG5_int tmp;
874
875 if ( !mesh->nt ) {
876 fprintf(stderr," ## Error: %s: You must set the number of elements with the",
877 __func__);
878 fprintf(stderr," MMG2D_Set_meshSize function before setting elements in mesh\n");
879 return 0;
880 }
881
882 if ( pos > mesh->ntmax ) {
883 fprintf(stderr," ## Error: %s: unable to allocate a new element.\n",
884 __func__);
885 fprintf(stderr," max number of element: %" MMG5_PRId "\n",mesh->ntmax);
887 return 0;
888 }
889
890 if ( pos > mesh->nt ) {
891 fprintf(stderr,"\n ## Error: %s: attempt to set new triangle at position %" MMG5_PRId ".",
892 __func__,pos);
893 fprintf(stderr," Overflow of the given number of triangle: %" MMG5_PRId "\n",mesh->nt);
894 fprintf(stderr," ## Check the mesh size, its compactness or the position");
895 fprintf(stderr," of the triangle.\n");
896 return 0;
897 }
898
899 pt = &mesh->tria[pos];
900 pt->v[0] = v0;
901 pt->v[1] = v1;
902 pt->v[2] = v2;
903 pt->ref = ref;
904
905 mesh->point[pt->v[0]].tag &= ~MG_NUL;
906 mesh->point[pt->v[1]].tag &= ~MG_NUL;
907 mesh->point[pt->v[2]].tag &= ~MG_NUL;
908
909 for(i=0 ; i<3 ; i++)
910 pt->edg[i] = 0;
911
912 vol = MMG2D_quickarea(mesh->point[pt->v[0]].c,mesh->point[pt->v[1]].c,
913 mesh->point[pt->v[2]].c);
914
915 if ( vol == 0.0 ) {
916 fprintf(stderr,"\n ## Error: %s: triangle %" MMG5_PRId " has null area.\n",
917 __func__,pos);
918 for ( ip=0; ip<3; ip++ ) {
919 ppt = &mesh->point[pt->v[ip]];
920 for ( j=0; j<3; j++ ) {
921 if ( fabs(ppt->c[j])>0. ) {
922 fprintf(stderr," Check that you don't have a sliver triangle.\n");
923 return 0;
924 }
925 }
926 }
927 }
928 else if(vol < 0) {
929 tmp = pt->v[2];
930 pt->v[2] = pt->v[1];
931 pt->v[1] = tmp;
932 /* mesh->xt temporary used to count reoriented tria */
933 mesh->xt++;
934 }
935 if ( mesh->info.ddebug && (mesh->nt == pos) && mesh->xt > 0 ) {
936 fprintf(stderr,"\n ## Warning: %s: %" MMG5_PRId " triangles reoriented\n",
937 __func__,mesh->xt);
938 mesh->xt = 0;
939 }
940
941 return 1;
942}
943
945 MMG5_pTria pt;
946
947 assert ( k <= mesh->nt );
948 pt = &mesh->tria[k];
949
950 pt->tag[0] |= MG_REQ;
951 pt->tag[1] |= MG_REQ;
952 pt->tag[2] |= MG_REQ;
953
954 return 1;
955}
956
958 MMG5_pTria pt;
959
960 assert ( k <= mesh->nt );
961 pt = &mesh->tria[k];
962
963 pt->tag[0] &= ~MG_REQ;
964 pt->tag[1] &= ~MG_REQ;
965 pt->tag[2] &= ~MG_REQ;
966
967 return 1;
968}
969
970int MMG2D_Get_triangle(MMG5_pMesh mesh, MMG5_int* v0, MMG5_int* v1, MMG5_int* v2, MMG5_int* ref
971 ,int* isRequired) {
972 MMG5_pTria ptt;
973
974 if ( mesh->nti == mesh->nt ) {
975 mesh->nti = 0;
976 if ( mesh->info.ddebug ) {
977 fprintf(stderr,"\n ## Warning: %s: reset the internal counter of"
978 " triangles.\n",__func__);
979 fprintf(stderr," You must pass here exactly one time (the first time ");
980 fprintf(stderr,"you call the MMG2D_Get_triangle function).\n");
981 fprintf(stderr," If not, the number of call of this function");
982 fprintf(stderr," exceed the number of triangles: %" MMG5_PRId "\n ",mesh->nt);
983 }
984 }
985
986 mesh->nti++;
987
988 if ( mesh->nti > mesh->nt ) {
989 fprintf(stderr,"\n ## Error: %s: unable to get triangle.\n",
990 __func__);
991 fprintf(stderr," The number of call of MMG2D_Get_triangle function");
992 fprintf(stderr," can not exceed the number of triangles: %" MMG5_PRId "\n ",mesh->nt);
993 return 0;
994 }
995
996 ptt = &mesh->tria[mesh->nti];
997
998 *v0 = ptt->v[0];
999 *v1 = ptt->v[1];
1000 *v2 = ptt->v[2];
1001 if ( ref != NULL )
1002 *ref = ptt->ref;
1003
1004 if ( isRequired != NULL ) {
1005 if ( (ptt->tag[0] & MG_REQ) && (ptt->tag[1] & MG_REQ) &&
1006 (ptt->tag[2] & MG_REQ) )
1007 *isRequired = 1;
1008 else
1009 *isRequired = 0;
1010 }
1011
1012 return 1;
1013}
1014
1015int MMG2D_Set_triangles(MMG5_pMesh mesh, MMG5_int *tria, MMG5_int *refs) {
1016 MMG5_pPoint ppt;
1017 MMG5_pTria ptt;
1018 double vol;
1019 int ip;
1020 MMG5_int j,i,tmp;
1021
1022 mesh->xt = 0;
1023 for (i=1;i<=mesh->nt;i++)
1024 {
1025 j = (i-1)*3;
1026 ptt = &mesh->tria[i];
1027 ptt->v[0] = tria[j] ;
1028 ptt->v[1] = tria[j+2];
1029 ptt->v[2] = tria[j+1];
1030 if ( refs != NULL )
1031 ptt->ref = refs[i-1];
1032
1033 mesh->point[ptt->v[0]].tag &= ~MG_NUL;
1034 mesh->point[ptt->v[1]].tag &= ~MG_NUL;
1035 mesh->point[ptt->v[2]].tag &= ~MG_NUL;
1036
1037 int ii;
1038 for( ii=0 ; ii<3 ; ii++)
1039 ptt->edg[ii] = 0;
1040
1041 vol = MMG2D_quickarea(mesh->point[ptt->v[0]].c,mesh->point[ptt->v[1]].c,
1042 mesh->point[ptt->v[2]].c);
1043
1044 if ( vol == 0.0 ) {
1045 fprintf(stderr,"\n ## Error: %s: triangle %" MMG5_PRId " has null area.\n",
1046 __func__,i);
1047 for ( ip=0; ip<3; ip++ ) {
1048 ppt = &mesh->point[ptt->v[ip]];
1049 int jj;
1050 for ( jj=0; jj<3; jj++ ) {
1051 if ( fabs(ppt->c[jj])>0. ) {
1052 fprintf(stderr," Check that you don't have a sliver triangle.\n");
1053 return 0;
1054 }
1055 }
1056 }
1057 }
1058 else if(vol < 0) {
1059 tmp = ptt->v[2];
1060 ptt->v[2] = ptt->v[1];
1061 ptt->v[1] = tmp;
1062 /* mesh->xt temporary used to count reoriented triangles */
1063 mesh->xt++;
1064 }
1065 if ( mesh->info.ddebug && mesh->xt > 0 ) {
1066 fprintf(stderr,"\n ## Warning: %s: %" MMG5_PRId " triangles reoriented\n",
1067 __func__,mesh->xt);
1068 }
1069 }
1070 return 1;
1071}
1072
1073int MMG2D_Get_triangles(MMG5_pMesh mesh, MMG5_int* tria, MMG5_int* refs,
1074 int* areRequired) {
1075 MMG5_pTria ptt;
1076 MMG5_int i, j;
1077
1078 for (i=1;i<=mesh->nt;i++)
1079 {
1080 j = (i-1)*3;
1081 ptt = &mesh->tria[i];
1082 tria[j] = ptt->v[0];
1083 tria[j+1] = ptt->v[1];
1084 tria[j+2] = ptt->v[2];
1085
1086 if ( refs!=NULL )
1087 refs[i-1] = ptt->ref ;
1088 if ( areRequired != NULL ) {
1089 if ( (ptt->tag[0] & MG_REQ) && (ptt->tag[1] & MG_REQ) &&
1090 (ptt->tag[2] & MG_REQ) )
1091 areRequired[i-1] = 1;
1092 else
1093 areRequired[i-1] = 0;
1094 }
1095 }
1096 return 1;
1097}
1098
1099int MMG2D_Set_quadrilateral(MMG5_pMesh mesh, MMG5_int v0, MMG5_int v1, MMG5_int v2, MMG5_int v3, MMG5_int ref, MMG5_int pos) {
1100 MMG5_pQuad pq;
1101
1102 if ( !mesh->nquad ) {
1103 fprintf(stderr,"\n ## Error: %s: You must set the number of quadrilaterals with the",
1104 __func__);
1105 fprintf(stderr," MMG2D_Set_meshSize function before setting elements in mesh\n");
1106 return 0;
1107 }
1108
1109 if ( pos > mesh->nquad ) {
1110 fprintf(stderr,"\n ## Error: %s: attempt to set new quad at position %" MMG5_PRId ".",
1111 __func__,pos);
1112 fprintf(stderr," Overflow of the given number of quads: %" MMG5_PRId "\n",mesh->nquad);
1113 fprintf(stderr,"\n ## Check the mesh size, its compactness or the position");
1114 fprintf(stderr," of the quad.\n");
1115 return 0;
1116 }
1117
1118 pq = &mesh->quadra[pos];
1119 pq->v[0] = v0;
1120 pq->v[1] = v1;
1121 pq->v[2] = v2;
1122 pq->v[3] = v3;
1123 pq->ref = ref;
1124
1125 mesh->point[pq->v[0]].tag &= ~MG_NUL;
1126 mesh->point[pq->v[1]].tag &= ~MG_NUL;
1127 mesh->point[pq->v[2]].tag &= ~MG_NUL;
1128 mesh->point[pq->v[3]].tag &= ~MG_NUL;
1129
1130 return 1;
1131}
1132
1133int MMG2D_Get_quadrilateral(MMG5_pMesh mesh, MMG5_int* v0, MMG5_int* v1, MMG5_int* v2, MMG5_int* v3,
1134 MMG5_int* ref, int* isRequired) {
1135 static MMG5_int nqi = 0;
1136
1137 if ( nqi == mesh->nquad ) {
1138 nqi = 0;
1139 if ( mesh->info.ddebug ) {
1140 fprintf(stderr,"\n ## Warning: %s: reset the internal counter of"
1141 " quadrilaterals.\n",__func__);
1142 fprintf(stderr," You must pass here exactly one time (the first time ");
1143 fprintf(stderr,"you call the MMG2D_Get_quadrilateral function).\n");
1144 fprintf(stderr," If not, the number of call of this function");
1145 fprintf(stderr," exceed the number of quadrilaterals: %" MMG5_PRId "\n ",mesh->nquad);
1146 }
1147 }
1148
1149 ++nqi;
1150
1151 if ( nqi > mesh->nquad ) {
1152 fprintf(stderr,"\n ## Error: %s: unable to get quadra.\n",__func__);
1153 fprintf(stderr," The number of call of MMG2D_Get_quadrilateral function");
1154 fprintf(stderr," can not exceed the number of quadra: %" MMG5_PRId "\n ",mesh->nquad);
1155 return 0;
1156 }
1157
1158 *v0 = mesh->quadra[nqi].v[0];
1159 *v1 = mesh->quadra[nqi].v[1];
1160 *v2 = mesh->quadra[nqi].v[2];
1161 *v3 = mesh->quadra[nqi].v[3];
1162
1163 if ( ref != NULL ) {
1164 *ref = mesh->quadra[nqi].ref;
1165 }
1166
1167 if ( isRequired != NULL ) {
1168 if ( ( mesh->quadra[nqi].tag[0] & MG_REQ) && ( mesh->quadra[nqi].tag[1] & MG_REQ) &&
1169 ( mesh->quadra[nqi].tag[2] & MG_REQ) && ( mesh->quadra[nqi].tag[3] & MG_REQ) ) {
1170 *isRequired = 1;
1171 }
1172 else
1173 *isRequired = 0;
1174 }
1175
1176 return 1;
1177}
1178
1179int MMG2D_Set_quadrilaterals(MMG5_pMesh mesh, MMG5_int *quadra, MMG5_int *refs) {
1180 MMG5_pQuad pq;
1181 MMG5_int i,j;
1182
1183 for (i=1;i<=mesh->nquad;i++)
1184 {
1185 j = (i-1)*4;
1186 pq = &mesh->quadra[i];
1187 pq->v[0] = quadra[j];
1188 pq->v[1] = quadra[j+1];
1189 pq->v[2] = quadra[j+2];
1190 pq->v[3] = quadra[j+3];
1191
1192 if ( refs != NULL )
1193 pq->ref = refs[i-1];
1194
1195 mesh->point[pq->v[0]].tag &= ~MG_NUL;
1196 mesh->point[pq->v[1]].tag &= ~MG_NUL;
1197 mesh->point[pq->v[2]].tag &= ~MG_NUL;
1198 mesh->point[pq->v[3]].tag &= ~MG_NUL;
1199 }
1200
1201 return 1;
1202}
1203
1204int MMG2D_Get_quadrilaterals(MMG5_pMesh mesh, MMG5_int *quadra, MMG5_int *refs, int * areRequired) {
1205 MMG5_pQuad pq;
1206 MMG5_int i, j;
1207
1208 for (i=1;i<=mesh->nquad;i++)
1209 {
1210 j = (i-1)*4;
1211 pq = &mesh->quadra[i];
1212 quadra[j] = pq->v[0];
1213 quadra[j+1] = pq->v[1];
1214 quadra[j+2] = pq->v[2];
1215 quadra[j+3] = pq->v[3];
1216
1217 if ( refs!=NULL )
1218 refs[i-1] = pq->ref ;
1219
1220 if ( areRequired != NULL ) {
1221 if ( (pq->tag[0] & MG_REQ) && (pq->tag[1] & MG_REQ) &&
1222 (pq->tag[2] & MG_REQ) && (pq->tag[3] & MG_REQ) ) {
1223 areRequired[i-1] = 1;
1224 }
1225 else
1226 areRequired[i-1] = 0;
1227 }
1228
1229 }
1230 return 1;
1231}
1232
1233
1234int MMG2D_Set_edge(MMG5_pMesh mesh, MMG5_int v0, MMG5_int v1, MMG5_int ref, MMG5_int pos) {
1235 MMG5_pEdge pt;
1236
1237 if ( !mesh->na ) {
1238 fprintf(stderr,"\n ## Error: %s: you must set the number of elements"
1239 " with the",__func__);
1240 fprintf(stderr," MMG2D_Set_meshSize function before setting elements in mesh\n");
1241 return 0;
1242 }
1243
1244 if ( pos > mesh->na ) {
1245 fprintf(stderr,"\n ## Error: %s: attempt to set new edge at position %" MMG5_PRId ".",
1246 __func__,pos);
1247 fprintf(stderr," Overflow of the given number of edge: %" MMG5_PRId "\n",mesh->na);
1248 fprintf(stderr," ## Check the mesh size, its compactness or the position");
1249 fprintf(stderr," of the edge.\n");
1250 return 0;
1251 }
1252
1253 pt = &mesh->edge[pos];
1254 pt->a = v0;
1255 pt->b = v1;
1256 pt->ref = ref;
1257 pt->tag &= MG_REF + MG_BDY;
1258
1259 mesh->point[pt->a].tag &= ~MG_NUL;
1260 mesh->point[pt->b].tag &= ~MG_NUL;
1261
1262 return 1;
1263}
1264
1266 MMG5_pEdge ped;
1267
1268 assert ( k <= mesh->na );
1269
1270 ped = &mesh->edge[k];
1271
1272 ped->tag |= MG_REQ;
1273
1274 return 1;
1275}
1276
1278 MMG5_pEdge ped;
1279
1280 assert ( k <= mesh->na );
1281
1282 ped = &mesh->edge[k];
1283
1284 ped->tag &= ~MG_REQ;
1285
1286 return 1;
1287}
1288
1290 MMG5_pPoint ppt;
1291 MMG5_pEdge ped;
1292
1293 assert ( k <= mesh->na );
1294
1295 ped = &mesh->edge[k];
1296
1297 ped->tag |= MG_PARBDY;
1298
1299 ppt = &mesh->point[ped->a];
1300 ppt->tag |= MG_PARBDY;
1301 ppt = &mesh->point[ped->b];
1302 ppt->tag |= MG_PARBDY;
1303
1304 return 1;
1305}
1306
1307int MMG2D_Get_edge(MMG5_pMesh mesh, MMG5_int* e0, MMG5_int* e1, MMG5_int* ref
1308 ,int* isRidge, int* isRequired) {
1309 MMG5_pEdge ped;
1310
1311 if ( mesh->nai == mesh->na ) {
1312 mesh->nai = 0;
1313 if ( mesh->info.ddebug ) {
1314 fprintf(stderr,"\n ## Warning: %s: reset the internal counter of edges.\n",
1315 __func__);
1316 fprintf(stderr," You must pass here exactly one time (the first time ");
1317 fprintf(stderr,"you call the MMG2D_Get_edge function).\n");
1318 fprintf(stderr," If not, the number of call of this function");
1319 fprintf(stderr," exceed the number of edges.\n ");
1320 fprintf(stderr," Please, call the MMG2D_Get_meshSize function to get"
1321 " this number.\n ");
1322 }
1323 }
1324
1325 mesh->nai++;
1326
1327 if ( mesh->nai > mesh->na ) {
1328 fprintf(stderr,"\n ## Error: %s: unable to get edge.\n",__func__);
1329 fprintf(stderr," The number of call of MMG2D_Get_edge function");
1330 fprintf(stderr," can not exceed the number of edges: %" MMG5_PRId "\n ",mesh->na);
1331 return 0;
1332 }
1333
1334 ped = &mesh->edge[mesh->nai];
1335
1336 while ( !ped->a && ++mesh->nai <= mesh->na ) {
1337 ped = &mesh->edge[mesh->nai];
1338 }
1339
1340 *e0 = ped->a;
1341 *e1 = ped->b;
1342
1343 if ( ref != NULL )
1344 *ref = mesh->edge[mesh->nai].ref;
1345
1346 if ( isRidge != NULL ) {
1347 if ( mesh->edge[mesh->nai].tag & MG_GEO )
1348 *isRidge = 1;
1349 else
1350 *isRidge = 0;
1351 }
1352
1353 if ( isRequired != NULL ) {
1354 if ( mesh->edge[mesh->nai].tag & MG_REQ )
1355 *isRequired = 1;
1356 else
1357 *isRequired = 0;
1358 }
1359
1360 return 1;
1361}
1362
1363int MMG2D_Set_edges(MMG5_pMesh mesh, MMG5_int *edges, MMG5_int *refs) {
1364 MMG5_int i,j;
1365
1366 for (i=1;i<=mesh->na;i++)
1367 {
1368 j = (i-1)*2;
1369
1370 mesh->edge[i].a = edges[j];
1371 mesh->edge[i].b = edges[j+1];
1372 if ( refs != NULL )
1373 mesh->edge[i].ref = refs[i-1];
1374 mesh->edge[i].tag &= MG_REF+MG_BDY;
1375
1376 mesh->point[mesh->edge[i].a].tag &= ~MG_NUL;
1377 mesh->point[mesh->edge[i].b].tag &= ~MG_NUL;
1378 }
1379
1380 return 1;
1381}
1382
1383int MMG2D_Get_edges(MMG5_pMesh mesh, MMG5_int* edges,MMG5_int *refs,int* areRidges,int* areRequired) {
1384 MMG5_int i,j;
1385
1386 for (i=1;i<=mesh->na;i++)
1387 {
1388 j = (i-1)*2;
1389 edges[j] = mesh->edge[i].a;
1390 edges[j+1] = mesh->edge[i].b;
1391
1392 if ( refs!=NULL )
1393 refs[i-1] = mesh->edge[i].ref;
1394
1395 if ( areRidges != NULL ) {
1396 if ( mesh->edge[i].tag & MG_GEO )
1397 areRidges[i-1] = 1;
1398 else
1399 areRidges[i-1] = 0;
1400 }
1401
1402 if ( areRequired != NULL ) {
1403 if ( mesh->edge[i].tag & MG_REQ )
1404 areRequired[i-1] = 1;
1405 else
1406 areRequired[i-1] = 0;
1407 }
1408 }
1409
1410 return 1;
1411}
1412
1414 double qual = 0.;
1415 MMG5_pTria pt;
1416
1417 if ( k < 1 || k > mesh->nt ) {
1418 fprintf(stderr,"\n ## Error: %s: unable to access to triangle %" MMG5_PRId ".\n",
1419 __func__,k);
1420 fprintf(stderr," Tria numbering goes from 1 to %" MMG5_PRId "\n",mesh->nt);
1421 return 0.;
1422 }
1423 pt = &mesh->tria[k];
1424 assert ( MG_EOK(pt) );
1425
1426 if ( (!met) || (!met->m) || met->size==1 ) {
1427 /* iso quality */
1428 qual = MMG2D_ALPHAD * MMG2D_caltri_iso(mesh,NULL,pt);
1429 }
1430 else {
1431 qual = MMG2D_ALPHAD * MMG2D_caltri_ani(mesh,met,pt);
1432 }
1433
1434 return qual;
1435}
1436
1437int MMG2D_Set_scalarSol(MMG5_pSol met, double s, MMG5_int pos) {
1438
1439 if ( !met->np ) {
1440 fprintf(stderr,"\n ## Error: %s: You must set the number of"
1441 " solution with the",__func__);
1442 fprintf(stderr," MMG2D_Set_solSize function before setting values");
1443 fprintf(stderr," in solution structure \n");
1444 return 0;
1445 }
1446
1447 if ( pos >= met->npmax ) {
1448 fprintf(stderr,"\n ## Error: %s: unable to set a new solution.\n",
1449 __func__);
1450 fprintf(stderr," max number of solutions: %" MMG5_PRId "\n",met->npmax);
1451 return 0;
1452 }
1453
1454 if ( pos > met->np ) {
1455 fprintf(stderr,"\n ## Error: %s: attempt to set new solution"
1456 " at position %" MMG5_PRId ".",__func__,pos);
1457 fprintf(stderr," Overflow of the given number of solutions: %" MMG5_PRId "\n",met->np);
1458 fprintf(stderr," ## Check the solution size, its compactness or the position");
1459 fprintf(stderr," of the solution.\n");
1460 return 0;
1461 }
1462
1463 met->m[pos] = s;
1464 return 1;
1465}
1466
1468{
1469 int ddebug = 0;
1470
1471 if ( met->npi == met->np ) {
1472 met->npi = 0;
1473 if ( ddebug ) {
1474 fprintf(stderr,"\n ## Warning: %s: reset the internal counter"
1475 " of points.\n",__func__);
1476 fprintf(stderr," You must pass here exactly one time (the first time ");
1477 fprintf(stderr,"you call the MMG2D_Get_scalarSol function).\n");
1478 fprintf(stderr," If not, the number of call of this function");
1479 fprintf(stderr," exceed the number of points: %" MMG5_PRId "\n ",met->np);
1480 }
1481 }
1482
1483 met->npi++;
1484
1485 if ( met->npi > met->np ) {
1486 fprintf(stderr,"\n ## Error: %s: unable to get solution.\n",
1487 __func__);
1488 fprintf(stderr," The number of call of MMG2D_Get_scalarSol function");
1489 fprintf(stderr," can not exceed the number of points: %" MMG5_PRId "\n ",met->np);
1490 return 0;
1491 }
1492
1493 *s = met->m[met->npi];
1494
1495 return 1;
1496}
1497
1498int MMG2D_Set_scalarSols(MMG5_pSol met, double *s) {
1499 MMG5_int k;
1500
1501 if ( !met->np ) {
1502 fprintf(stderr,"\n ## Error: %s: you must set the number of"
1503 " solution with the",__func__);
1504 fprintf(stderr," MMG2D_Set_solSize function before setting values");
1505 fprintf(stderr," in solution structure \n");
1506 return 0;
1507 }
1508
1509 for ( k=0; k<met->np; ++k )
1510 met->m[k+1] = s[k];
1511
1512 return 1;
1513}
1514
1515int MMG2D_Get_scalarSols(MMG5_pSol met, double* s) {
1516 MMG5_int k;
1517
1518 for ( k=0; k<met->np; ++k )
1519 s[k] = met->m[k+1];
1520
1521 return 1;
1522}
1523
1524int MMG2D_Set_vectorSol(MMG5_pSol met, double vx,double vy, MMG5_int pos) {
1525 MMG5_int isol;
1526
1527 if ( !met->np ) {
1528 fprintf(stderr,"\n ## Error: %s: you must set the number of"
1529 " solution with the",__func__);
1530 fprintf(stderr," MMG2D_Set_solSize function before setting values");
1531 fprintf(stderr," in solution structure \n");
1532 return 0;
1533 }
1534 if ( pos < 1 ) {
1535 fprintf(stderr,"\n ## Error: %s: unable to set a new solution.\n",
1536 __func__);
1537 fprintf(stderr," Minimal index of the solution position must be 1.\n");
1538 return 0;
1539 }
1540 if ( pos >= met->npmax ) {
1541 fprintf(stderr,"\n ## Error: %s: unable to set a new solution.\n",
1542 __func__);
1543 fprintf(stderr," max number of solutions: %" MMG5_PRId "\n",met->npmax);
1544 return 0;
1545 }
1546
1547 if ( pos > met->np ) {
1548 fprintf(stderr,"\n ## Error: %s: attempt to set new solution"
1549 " at position %" MMG5_PRId ".",__func__,pos);
1550 fprintf(stderr," Overflow of the given number of solutions: %" MMG5_PRId "\n",met->np);
1551 fprintf(stderr,"\n ## Check the solution size, its compactness or the position");
1552 fprintf(stderr," of the solution.\n");
1553 return 0;
1554 }
1555
1556 isol = (pos-1) * met->size + 1;
1557
1558 met->m[isol] = vx;
1559 met->m[isol+1] = vy;
1560
1561 return 1;
1562}
1563
1564
1565int MMG2D_Get_vectorSol(MMG5_pSol met, double* vx, double* vy) {
1566
1567 int ddebug = 0;
1568
1569 if ( met->npi == met->np ) {
1570 met->npi = 0;
1571 if ( ddebug ) {
1572 fprintf(stderr,"\n ## Warning: %s: reset the internal counter of points.\n",
1573 __func__);
1574 fprintf(stderr," You must pass here exactly one time (the first time ");
1575 fprintf(stderr,"you call the MMG2D_Get_vectorSol function).\n");
1576 fprintf(stderr," If not, the number of call of this function");
1577 fprintf(stderr," exceed the number of points: %" MMG5_PRId "\n ",met->np);
1578 }
1579 }
1580
1581 met->npi++;
1582
1583 if ( met->npi > met->np ) {
1584 fprintf(stderr,"\n ## Error: %s: unable to get solution.\n",__func__);
1585 fprintf(stderr," The number of call of MMG2D_Get_vectorSol function");
1586 fprintf(stderr," can not exceed the number of points: %" MMG5_PRId "\n ",met->np);
1587 return 0;
1588 }
1589
1590 *vx = met->m[met->size*(met->npi-1)+1];
1591 *vy = met->m[met->size*(met->npi-1)+2];
1592
1593 return 1;
1594}
1595
1596int MMG2D_Set_vectorSols(MMG5_pSol met, double *sols) {
1597 double *m;
1598 MMG5_int k,j;
1599
1600 if ( !met->np ) {
1601 fprintf(stderr,"\n ## Error: %s: you must set the number of"
1602 " solution with the",__func__);
1603 fprintf(stderr," MMG2D_Set_solSize function before setting values");
1604 fprintf(stderr," in solution structure \n");
1605 return 0;
1606 }
1607
1608 for ( k=0; k<met->np; ++k ) {
1609 j = 2*k;
1610 m = &met->m[j];
1611 m[1] = sols[j];
1612 m[2] = sols[j+1];
1613 }
1614
1615 return 1;
1616}
1617
1618int MMG2D_Get_vectorSols(MMG5_pSol met, double* sols) {
1619 double *m;
1620 MMG5_int k, j;
1621
1622 for ( k=0; k<met->np; ++k ) {
1623 j = 2*k;
1624 m = &met->m[j];
1625
1626 sols[j] = m[1];
1627 sols[j+1] = m[2];
1628 }
1629
1630 return 1;
1631}
1632
1633
1634int MMG2D_Set_tensorSol(MMG5_pSol met, double m11, double m12, double m22,
1635 MMG5_int pos) {
1636 MMG5_int isol;
1637
1638 if ( !met->np ) {
1639 fprintf(stderr,"\n ## Error: %s: you must set the number of"
1640 " solution with the",__func__);
1641 fprintf(stderr," MMG2D_Set_solSize function before setting values");
1642 fprintf(stderr," in solution structure \n");
1643 return 0;
1644 }
1645
1646 if ( pos >= met->npmax ) {
1647 fprintf(stderr,"\n ## Error: %s: unable to set a new solution.\n",
1648 __func__);
1649 fprintf(stderr," max number of solutions: %" MMG5_PRId "\n",met->npmax);
1650 return 0;
1651 }
1652
1653 if ( pos > met->np ) {
1654 fprintf(stderr,"\n ## Error: %s: attempt to set new solution "
1655 "at position %" MMG5_PRId ".",__func__,pos);
1656 fprintf(stderr," Overflow of the given number of solutions: %" MMG5_PRId "\n",met->np);
1657 fprintf(stderr," ## Check the solution size, its compactness or the position");
1658 fprintf(stderr," of the solution.\n");
1659 return 0;
1660 }
1661 isol = pos * met->size;
1662 met->m[isol ] = m11;
1663 met->m[isol + 1] = m12;
1664 met->m[isol + 2] = m22;
1665 return 1;
1666}
1667
1668int MMG2D_Get_tensorSol(MMG5_pSol met, double *m11,double *m12,double *m22)
1669{
1670 int ddebug = 0;
1671 MMG5_int isol;
1672
1673 if ( met->npi == met->np ) {
1674 met->npi = 0;
1675 if ( ddebug ) {
1676 fprintf(stderr,"\n ## Warning: %s: reset the internal counter of points.\n",
1677 __func__);
1678 fprintf(stderr," You must pass here exactly one time (the first time ");
1679 fprintf(stderr,"you call the MMG2D_Get_tensorSol function).\n");
1680 fprintf(stderr," If not, the number of call of this function");
1681 fprintf(stderr," exceed the number of points: %" MMG5_PRId "\n ",met->np);
1682 }
1683 }
1684
1685 met->npi++;
1686
1687 if ( met->npi > met->np ) {
1688 fprintf(stderr,"\n ## Error: %s: unable to get solution.\n",__func__);
1689 fprintf(stderr," The number of call of MMG2D_Get_tensorSol function");
1690 fprintf(stderr," can not exceed the number of points: %" MMG5_PRId "\n ",met->np);
1691 return 0;
1692 }
1693
1694 isol = met->size*(met->npi);
1695 *m11 = met->m[isol ];
1696 *m12 = met->m[isol+1];
1697 *m22 = met->m[isol+2];
1698
1699 return 1;
1700}
1701
1702int MMG2D_Set_tensorSols(MMG5_pSol met, double *sols) {
1703 double *m;
1704 MMG5_int k,j;
1705
1706 if ( !met->np ) {
1707 fprintf(stderr,"\n ## Error: %s: You must set the number"
1708 " of solution with the",__func__);
1709 fprintf(stderr," MMG2D_Set_solSize function before setting values");
1710 fprintf(stderr," in solution structure \n");
1711 return 0;
1712 }
1713
1714 for ( k=0; k<met->np; ++k ) {
1715 j = 3*k;
1716 m = &met->m[j+3];
1717
1718 m[0] = sols[j];
1719 m[1] = sols[j+1];
1720 m[2] = sols[j+2];
1721 }
1722 return 1;
1723}
1724
1725int MMG2D_Get_tensorSols(MMG5_pSol met, double *sols) {
1726 double *m;
1727 MMG5_int k,j;
1728
1729 for ( k=0; k<met->np; ++k ) {
1730 j = 3*k;
1731 m = &met->m[j+3];
1732
1733 sols[j] = m[0];
1734 sols[j+1] = m[1];
1735 sols[j+2] = m[2];
1736 }
1737
1738 return 1;
1739}
1740
1742 MMG5_pSol psl;
1743
1744 /* Warning: users give indices from 1 to nsols */
1745 psl = sol + (i-1);
1746
1747 switch ( psl->type ) {
1748 case MMG5_Scalar:
1749 return MMG2D_Set_scalarSols(psl,s);
1750 break;
1751
1752 case MMG5_Vector:
1753 MMG2D_Set_vectorSols(psl,s);
1754 break;
1755
1756 case MMG5_Tensor:
1757 MMG2D_Set_tensorSols(psl,s);
1758 break;
1759
1760 default:
1761 fprintf(stderr,"\n ## Error: %s: unexpected type of solution: %s.\n",
1762 __func__,MMG5_Get_typeName(psl->type));
1763 return 0;
1764 }
1765
1766 return 1;
1767}
1768
1770 MMG5_pSol psl;
1771
1772 /* Warning: users give indices from 1 to nsols */
1773 psl = sol + (i-1);
1774
1775 switch ( psl->type ) {
1776 case MMG5_Scalar:
1777 return MMG2D_Get_scalarSols(psl,s);
1778 break;
1779
1780 case MMG5_Vector:
1781 MMG2D_Get_vectorSols(psl,s);
1782 break;
1783
1784 case MMG5_Tensor:
1785 MMG2D_Get_tensorSols(psl,s);
1786 break;
1787
1788 default:
1789 fprintf(stderr,"\n ## Error: %s: unexpected type of solution: %s\n",
1790 __func__,MMG5_Get_typeName(psl->type));
1791 return 0;
1792 }
1793
1794 return 1;
1795}
1796
1797int MMG2D_Set_ithSol_inSolsAtVertices(MMG5_pSol sol,int i, double* s,MMG5_int pos) {
1798 MMG5_pSol psl;
1799
1800 /* Warning: users give indices from 1 to nsols */
1801 psl = sol + (i-1);
1802
1803 switch ( psl->type ) {
1804 case MMG5_Scalar:
1805 return MMG2D_Set_scalarSol(psl,s[0],pos);
1806 break;
1807
1808 case MMG5_Vector:
1809 MMG2D_Set_vectorSol(psl,s[0],s[1],pos);
1810 break;
1811
1812 case MMG5_Tensor:
1813 MMG2D_Set_tensorSol(psl,s[0],s[1],s[2],pos);
1814 break;
1815
1816 default:
1817 fprintf(stderr,"\n ## Error: %s: unexpected type of solution: %s.\n",
1818 __func__,MMG5_Get_typeName(psl->type));
1819 return 0;
1820 }
1821 return 1;
1822}
1823
1824int MMG2D_Get_ithSol_inSolsAtVertices(MMG5_pSol sol,int i, double *s,MMG5_int pos) {
1825 MMG5_pSol psl;
1826
1827 /* Warning: users give indices from 1 to nsols */
1828 psl = sol + (i-1);
1829
1830 psl->npi = pos-1;
1831
1832 switch ( psl->type ) {
1833 case MMG5_Scalar:
1834 return MMG2D_Get_scalarSol(psl,&s[0]);
1835 break;
1836
1837 case MMG5_Vector:
1838 MMG2D_Get_vectorSol(psl,&s[0],&s[1]);
1839 break;
1840
1841 case MMG5_Tensor:
1842 MMG2D_Get_tensorSol(psl,&s[0],&s[1],&s[2]);
1843 break;
1844
1845 default:
1846 fprintf(stderr,"\n ## Error: %s: unexpected type of solution: %s\n",
1847 __func__,MMG5_Get_typeName(psl->type));
1848 return 0;
1849 }
1850
1851 return 1;
1852}
1853
1855
1856 if ( (mesh->npi != mesh->np) || (mesh->nti != mesh->nt) ) {
1857 fprintf(stderr,"\n ## Error: %s: if you don't use the MMG2D_loadMesh function,",
1858 __func__);
1859 fprintf(stderr," you must call the MMG2D_Set_meshSize function to have a");
1860 fprintf(stderr," valid mesh.\n");
1861 fprintf(stderr," Missing datas.\n");
1862 return 0;
1863 }
1864
1865 if ( met->npi != met->np ) {
1866 fprintf(stderr,"\n ## Error: %s: if you don't use the MMG2D_loadMet function,",
1867 __func__);
1868 fprintf(stderr," you must call the MMG2D_Set_solSize function to have a");
1869 fprintf(stderr," valid solution.\n");
1870 fprintf(stderr," Missing datas.\n");
1871 return 0;
1872 }
1873
1874 /* Check mesh data */
1875 if ( mesh->info.ddebug ) {
1876 if ( (!mesh->np) || (!mesh->point) ||
1877 (!mesh->nt) ) {
1878 fprintf(stderr," ** MISSING DATA.\n");
1879 fprintf(stderr," Check that your mesh contains points.\n");
1880 fprintf(stderr," Exit program.\n");
1881 return 0;
1882 }
1883 }
1884
1885 if ( mesh->dim != 2 ) {
1886 fprintf(stderr," ** 2 DIMENSIONAL MESH NEEDED. Exit program.\n");
1887 return 0;
1888 }
1889 if ( met->dim != 2 ) {
1890 fprintf(stderr," ** WRONG DIMENSION FOR METRIC. Exit program.\n");
1891 return 0;
1892 }
1893 if ( !mesh->ver ) mesh->ver = 2;
1894 if ( !met ->ver ) met ->ver = 2;
1895
1896 return 1;
1897}
1898
1903
1904
1905int MMG2D_Free_all(const int starter,...)
1906{
1907 va_list argptr;
1908 int ier;
1909
1911
1913
1914 va_end(argptr);
1915
1916 return ier;
1917}
1918
1920{
1921 int ier;
1922
1923 va_list argptr;
1924
1926
1928
1929 va_end(argptr);
1930
1931 return ier;
1932}
1933
1934int MMG2D_Free_names(const int starter,...)
1935{
1936 va_list argptr;
1937 int ier;
1938
1940
1942
1943 va_end(argptr);
1944
1945 return ier;
1946}
const char * MMG5_Get_typeName(enum MMG5_type typ)
int MMG5_Set_outputMeshName(MMG5_pMesh mesh, const char *meshout)
int MMG5_Set_inputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solin)
int MMG5_Free_allSols(MMG5_pMesh mesh, MMG5_pSol *sol)
int MMG5_Set_inputMeshName(MMG5_pMesh mesh, const char *meshin)
void MMG5_Init_parameters(MMG5_pMesh mesh)
const char * MMG5_Get_entitiesName(enum MMG5_entities ent)
void MMG5_Init_fileNames(MMG5_pMesh mesh, MMG5_pSol sol)
int MMG5_Set_outputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solout)
int MMG2D_Set_lsBaseReference(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int br)
int MMG2D_Set_scalarSol(MMG5_pSol met, double s, MMG5_int pos)
int MMG2D_Get_ithSols_inSolsAtVertices(MMG5_pSol sol, int i, double *s)
void MMG2D_Init_parameters(MMG5_pMesh mesh)
void MMG2D_Init_fileNames(MMG5_pMesh mesh, MMG5_pSol sol)
int MMG2D_Set_quadrilateral(MMG5_pMesh mesh, MMG5_int v0, MMG5_int v1, MMG5_int v2, MMG5_int v3, MMG5_int ref, MMG5_int pos)
int MMG2D_Unset_requiredEdge(MMG5_pMesh mesh, MMG5_int k)
int MMG2D_Set_vertex(MMG5_pMesh mesh, double c0, double c1, MMG5_int ref, MMG5_int pos)
int MMG2D_Set_scalarSols(MMG5_pSol met, double *s)
int MMG2D_Init_mesh(const int starter,...)
int MMG2D_Get_vectorSols(MMG5_pSol met, double *sols)
int MMG2D_Set_iparameter(MMG5_pMesh mesh, MMG5_pSol sol, int iparam, MMG5_int val)
int MMG2D_Set_meshSize(MMG5_pMesh mesh, MMG5_int np, MMG5_int nt, MMG5_int nquad, MMG5_int na)
int MMG2D_Set_requiredEdge(MMG5_pMesh mesh, MMG5_int k)
int MMG2D_Get_quadrilateral(MMG5_pMesh mesh, MMG5_int *v0, MMG5_int *v1, MMG5_int *v2, MMG5_int *v3, MMG5_int *ref, int *isRequired)
int MMG2D_Set_tensorSol(MMG5_pSol met, double m11, double m12, double m22, MMG5_int pos)
int MMG2D_Get_vectorSol(MMG5_pSol met, double *vx, double *vy)
int MMG2D_Set_edge(MMG5_pMesh mesh, MMG5_int v0, MMG5_int v1, MMG5_int ref, MMG5_int pos)
int MMG2D_Set_ithSols_inSolsAtVertices(MMG5_pSol sol, int i, double *s)
int MMG2D_Unset_requiredTriangle(MMG5_pMesh mesh, MMG5_int k)
int MMG2D_Set_triangles(MMG5_pMesh mesh, MMG5_int *tria, MMG5_int *refs)
int MMG2D_GetByIdx_vertex(MMG5_pMesh mesh, double *c0, double *c1, MMG5_int *ref, int *isCorner, int *isRequired, MMG5_int idx)
int MMG2D_Free_allSols(MMG5_pMesh mesh, MMG5_pSol *sol)
int MMG2D_Set_ithSol_inSolsAtVertices(MMG5_pSol sol, int i, double *s, MMG5_int pos)
int MMG2D_Get_scalarSol(MMG5_pSol met, double *s)
int MMG2D_Get_vertices(MMG5_pMesh mesh, double *vertices, MMG5_int *refs, int *areCorners, int *areRequired)
int MMG2D_Unset_corner(MMG5_pMesh mesh, MMG5_int k)
int MMG2D_Get_solSize(MMG5_pMesh mesh, MMG5_pSol sol, int *typEntity, MMG5_int *np, int *typSol)
int MMG2D_Get_tensorSol(MMG5_pSol met, double *m11, double *m12, double *m22)
int MMG2D_Set_solSize(MMG5_pMesh mesh, MMG5_pSol sol, int typEntity, MMG5_int np, int typSol)
int MMG2D_Free_structures(const int starter,...)
int MMG2D_Chk_meshData(MMG5_pMesh mesh, MMG5_pSol met)
int MMG2D_Set_inputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solin)
int MMG2D_Set_tensorSols(MMG5_pSol met, double *sols)
int MMG2D_Set_inputMeshName(MMG5_pMesh mesh, const char *meshin)
int MMG2D_Get_vertex(MMG5_pMesh mesh, double *c0, double *c1, MMG5_int *ref, int *isCorner, int *isRequired)
int MMG2D_Get_scalarSols(MMG5_pSol met, double *s)
int MMG2D_Set_vectorSol(MMG5_pSol met, double vx, double vy, MMG5_int pos)
int MMG2D_Get_edge(MMG5_pMesh mesh, MMG5_int *e0, MMG5_int *e1, MMG5_int *ref, int *isRidge, int *isRequired)
double MMG2D_Get_triangleQuality(MMG5_pMesh mesh, MMG5_pSol met, MMG5_int k)
int MMG2D_Set_triangle(MMG5_pMesh mesh, MMG5_int v0, MMG5_int v1, MMG5_int v2, MMG5_int ref, MMG5_int pos)
int MMG2D_Set_solsAtVerticesSize(MMG5_pMesh mesh, MMG5_pSol *sol, int nsols, MMG5_int np, int *typSol)
int MMG2D_Get_tensorSols(MMG5_pSol met, double *sols)
int MMG2D_Get_meshSize(MMG5_pMesh mesh, MMG5_int *np, MMG5_int *nt, MMG5_int *nquad, MMG5_int *na)
int MMG2D_Get_edges(MMG5_pMesh mesh, MMG5_int *edges, MMG5_int *refs, int *areRidges, int *areRequired)
int MMG2D_Get_quadrilaterals(MMG5_pMesh mesh, MMG5_int *quadra, MMG5_int *refs, int *areRequired)
int MMG2D_Set_localParameter(MMG5_pMesh mesh, MMG5_pSol sol, int typ, MMG5_int ref, double hmin, double hmax, double hausd)
int MMG2D_Set_vertices(MMG5_pMesh mesh, double *vertices, MMG5_int *refs)
int MMG2D_Get_triangle(MMG5_pMesh mesh, MMG5_int *v0, MMG5_int *v1, MMG5_int *v2, MMG5_int *ref, int *isRequired)
int MMG2D_Set_vectorSols(MMG5_pSol met, double *sols)
int MMG2D_Set_corner(MMG5_pMesh mesh, MMG5_int k)
int MMG2D_Set_edges(MMG5_pMesh mesh, MMG5_int *edges, MMG5_int *refs)
int MMG2D_Set_requiredVertex(MMG5_pMesh mesh, MMG5_int k)
int MMG2D_Get_ithSol_inSolsAtVertices(MMG5_pSol sol, int i, double *s, MMG5_int pos)
int MMG2D_Set_parallelEdge(MMG5_pMesh mesh, MMG5_int k)
int MMG2D_Set_dparameter(MMG5_pMesh mesh, MMG5_pSol sol, int dparam, double val)
int MMG2D_Free_all(const int starter,...)
int MMG2D_Set_multiMat(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int ref, int split, MMG5_int rin, MMG5_int rout)
int MMG2D_Get_solsAtVerticesSize(MMG5_pMesh mesh, MMG5_pSol *sol, int *nsols, MMG5_int *np, int *typSol)
int MMG2D_Get_triangles(MMG5_pMesh mesh, MMG5_int *tria, MMG5_int *refs, int *areRequired)
int MMG2D_Set_quadrilaterals(MMG5_pMesh mesh, MMG5_int *quadra, MMG5_int *refs)
int MMG2D_Free_names(const int starter,...)
int MMG2D_Unset_requiredVertex(MMG5_pMesh mesh, MMG5_int k)
int MMG2D_Set_outputMeshName(MMG5_pMesh mesh, const char *meshout)
int MMG2D_Set_requiredTriangle(MMG5_pMesh mesh, MMG5_int k)
int MMG2D_Set_outputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solout)
int ier
tmp[*strlen0]
const int starter
MMG5_pMesh MMG5_pSol * sol
va_start(argptr, starter)
MMG5_pMesh char * meshin
MMG5_pMesh * mesh
va_end(argptr)
const int va_list argptr
API headers for the mmg2d library.
@ MMG2D_DPARAM_hgradreq
Definition libmmg2d.h:87
@ MMG2D_IPARAM_numsubdomain
Definition libmmg2d.h:75
@ MMG2D_IPARAM_iso
Definition libmmg2d.h:63
@ MMG2D_IPARAM_optim
Definition libmmg2d.h:68
@ MMG2D_IPARAM_numberOfLocalParam
Definition libmmg2d.h:76
@ MMG2D_IPARAM_nosurf
Definition libmmg2d.h:72
@ MMG2D_DPARAM_hgrad
Definition libmmg2d.h:86
@ MMG2D_DPARAM_hmin
Definition libmmg2d.h:82
@ MMG2D_IPARAM_mem
Definition libmmg2d.h:60
@ MMG2D_IPARAM_nofem
Definition libmmg2d.h:90
@ MMG2D_IPARAM_nosizreq
Definition libmmg2d.h:80
@ MMG2D_IPARAM_isoref
Definition libmmg2d.h:91
@ MMG2D_DPARAM_rmc
Definition libmmg2d.h:89
@ MMG2D_DPARAM_hausd
Definition libmmg2d.h:85
@ MMG2D_IPARAM_angle
Definition libmmg2d.h:62
@ MMG2D_DPARAM_hmax
Definition libmmg2d.h:83
@ MMG2D_IPARAM_isosurf
Definition libmmg2d.h:64
@ MMG2D_DPARAM_ls
Definition libmmg2d.h:88
@ MMG2D_IPARAM_noinsert
Definition libmmg2d.h:69
@ MMG2D_IPARAM_nreg
Definition libmmg2d.h:73
@ MMG2D_IPARAM_noswap
Definition libmmg2d.h:70
@ MMG2D_IPARAM_xreg
Definition libmmg2d.h:74
@ MMG2D_IPARAM_nomove
Definition libmmg2d.h:71
@ MMG2D_IPARAM_lag
Definition libmmg2d.h:66
@ MMG2D_IPARAM_opnbdy
Definition libmmg2d.h:65
@ MMG2D_IPARAM_verbose
Definition libmmg2d.h:59
@ MMG2D_IPARAM_3dMedit
Definition libmmg2d.h:67
@ MMG2D_IPARAM_anisosize
Definition libmmg2d.h:79
@ MMG2D_IPARAM_numberOfLSBaseReferences
Definition libmmg2d.h:77
@ MMG2D_IPARAM_numberOfMat
Definition libmmg2d.h:78
@ MMG2D_DPARAM_angleDetection
Definition libmmg2d.h:81
@ MMG2D_DPARAM_hsiz
Definition libmmg2d.h:84
@ MMG2D_IPARAM_debug
Definition libmmg2d.h:61
int MMG2D_setMeshSize_alloc(MMG5_pMesh)
Definition zaldy_2d.c:252
int MMG2D_memOption(MMG5_pMesh mesh)
Definition zaldy_2d.c:233
double MMG2D_caltri_ani(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_pTria)
Definition quality_2d.c:121
double MMG2D_caltri_iso(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_pTria)
Definition quality_2d.c:107
#define MMG2D_ALPHAD
int MMG2D_Free_all_var(va_list argptr)
int MMG2D_Free_structures_var(va_list argptr)
int MMG2D_Free_names_var(va_list argptr)
int MMG2D_Init_mesh_var(va_list argptr)
LIBMMG_CORE_EXPORT int MMG5_Set_multiMat(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int ref, int split, MMG5_int rin, MMG5_int rex)
Definition libtools.c:102
#define MMG5_VOLFRAC
API header for the common part of the MMG libraries.
LIBMMG_CORE_EXPORT int MMG5_Set_lsBaseReference(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int br)
Definition libtools.c:174
@ MMG5_Vector
@ MMG5_Tensor
@ MMG5_Scalar
@ MMG5_Notype
@ MMG5_Noentity
@ MMG5_Vertex
@ MMG5_Edg
@ MMG5_Triangle
#define MG_Tria
#define MG_Edge
#define MG_REQ
#define MG_GEO
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
double MMG2D_quickarea(double a[2], double b[2], double c[2])
Definition tools.c:971
#define MG_EOK(pt)
#define MG_NUL
#define MMG5_INCREASE_MEM_MESSAGE()
#define MG_PARBDY
#define MMG5_LAG
#define MG_MIN(a, b)
#define MG_MAX(a, b)
#define MMG5_ADD_MEM(mesh, size, message, law)
#define MMG5_LS
#define MMG5_OFF
#define MG_CRN
#define MMG5_NR
#define MG_BDY
#define MMG5_ANGEDG
#define MMG5_NOHGRAD
#define MG_REF
#define MMG5_FEM
#define MMG5_DEL_MEM(mesh, ptr)
Structure to store edges of a MMG mesh.
MMG5_int b
MMG5_int ref
int16_t tag
MMG5_int a
MMG5_hgeom * geom
int8_t iso
int8_t parTyp
int8_t ddebug
double hsiz
int8_t isosurf
int8_t sethmin
MMG5_int * br
uint8_t ani
uint8_t noswap
double rmc
double hmin
int8_t setfem
MMG5_pMat mat
double hgrad
uint8_t noinsert
MMG5_int isoref
double hmax
double ls
uint8_t nomove
int8_t lag
MMG5_int nsd
uint8_t nosurf
int8_t sethmax
uint8_t nosizreq
MMG5_pPar par
uint8_t optim
double dhd
double hgradreq
double hausd
int8_t xreg
int8_t nreg
To store user-defined references in the mesh (useful in LS mode)
MMG mesh structure.
MMG5_int ntmax
MMG5_pQuad quadra
MMG5_Info info
MMG5_int xt
MMG5_pPoint point
MMG5_int nquad
MMG5_int npmax
MMG5_pxPoint xpoint
MMG5_HGeom htab
MMG5_int namax
MMG5_int nt
MMG5_pTria tria
MMG5_int np
MMG5_pEdge edge
MMG5_int nti
MMG5_int npi
MMG5_pxTetra xtetra
MMG5_int nai
MMG5_int na
double hmin
double hmax
double hausd
MMG5_int ref
int8_t elt
Structure to store points of a MMG mesh.
int16_t tag
MMG5_int tmp
double c[3]
MMG5_int ref
MMG5_int flag
MMG5_int ref
int16_t tag[4]
MMG5_int v[4]
MMG5_int npi
MMG5_int npmax
double * m
MMG5_int np
MMG5_int edg[3]
int16_t tag[3]
MMG5_int ref
MMG5_int v[3]