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