Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
inout_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
24#include "libmmg2d.h"
25#include "libmmg2d_private.h"
26
27/* read mesh data */
29 FILE *inm;
30 MMG5_pPoint ppt;
31 MMG5_pEdge ped;
32 MMG5_pTria pt;
33 MMG5_pQuad pq1;
34 float fc;
35 long posnp,posnt,posncor,posned,posnq,posreq,posreqed,posntreq,posnqreq;
36 MMG5_int k,tmp,ncor,norient,nreq,ntreq,nreqed,nqreq,nref;
37 int bin,iswp;
38 double air,dtmp;
39 int bdim,binch,bpos;
40 MMG5_int ref,i;
41 char *ptr,*data;
42 char chaine[MMG5_FILESTR_LGTH],strskip[MMG5_FILESTR_LGTH];
43
44 posnp = posnt = posncor = posned = posnq = posreq = posreqed = posntreq = posnqreq = 0;
45 ncor = nreq = nreqed = ntreq = nqreq = 0;
46 bin = 0;
47 iswp = 0;
48 mesh->np = mesh->nt = mesh->na = mesh->xp = 0;
49 nref = 0;
50
51 MMG5_SAFE_CALLOC(data,strlen(filename)+7,char,return -1);
52 strcpy(data,filename);
53 ptr = strstr(data,".mesh");
54 if ( !ptr ) {
55 strcat(data,".meshb");
56 if (!(inm = fopen(data,"rb")) ) {
57 ptr = strstr(data,".mesh");
58 *ptr = '\0';
59 strcat(data,".mesh");
60 if (!(inm = fopen(data,"rb")) ) {
61 MMG5_SAFE_FREE(data);
62 return 0;
63 }
64 }
65 else bin = 1;
66 }
67 else {
68 ptr = strstr(data,".meshb");
69
70 if ( ptr ) bin = 1;
71
72 if( !(inm = fopen(data,"rb")) ) {
73 MMG5_SAFE_FREE(data);
74 return 0;
75 }
76 }
77 if ( mesh->info.imprim >= 0 ) {
78 fprintf(stdout," %%%% %s OPENED\n",data);
79 }
80 MMG5_SAFE_FREE(data);
81
82 if (!bin) {
83 strcpy(chaine,"D");
84 while(fscanf(inm,"%127s",&chaine[0])!=EOF && strncmp(chaine,"End",strlen("End")) ) {
85 if ( chaine[0] == '#' ) {
86 fgets(strskip,MMG5_FILESTR_LGTH,inm);
87 continue;
88 }
89
90 if(!strncmp(chaine,"MeshVersionFormatted",strlen("MeshVersionFormatted"))) {
91 MMG_FSCANF(inm,"%d",&mesh->ver);
92 continue;
93 }
94 else if(!strncmp(chaine,"Dimension",strlen("Dimension"))) {
95 MMG_FSCANF(inm,"%d",&mesh->dim);
96 if ( mesh->info.renum >= 2) {
97 if(mesh->dim!=3) {
98 fprintf(stdout,"WRONG USE OF 3dMedit option \n");
99 return -1;
100 }
101 mesh->dim = 2;
102 }
103 if(mesh->dim!=2) {
104 fprintf(stdout,"BAD DIMENSION : %d\n",mesh->dim);
105 return -1;
106 }
107 continue;
108 }
109 else if(!strncmp(chaine,"Vertices",strlen("Vertices"))) {
110 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->np);
111 posnp = ftell(inm);
112 continue;
113 }
114 else if(!strncmp(chaine,"Triangles",strlen("Triangles"))) {
115 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->nt);
116 posnt = ftell(inm);
117 continue;
118 }
119 else if(!strncmp(chaine,"Quadrilaterals",strlen("Quadrilaterals"))) {
120 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->nquad);
121 posnq = ftell(inm);
122 continue;
123 }
124 else if(!strncmp(chaine,"RequiredQuadrilaterals",strlen("RequiredQuadrilaterals"))) {
125 MMG_FSCANF(inm,"%" MMG5_PRId "",&nqreq);
126 posnqreq = ftell(inm);
127 continue;
128 }
129 else if(!strncmp(chaine,"Corners",strlen("Corners"))) {
130 MMG_FSCANF(inm,"%" MMG5_PRId "",&ncor);
131 posncor = ftell(inm);
132 continue;
133 }
134 else if(!strncmp(chaine,"RequiredVertices",strlen("RequiredVertices"))) {
135 MMG_FSCANF(inm,"%" MMG5_PRId "",&nreq);
136 posreq = ftell(inm);
137 continue;
138 }
139 else if(!strncmp(chaine,"Edges",strlen("Edges"))) {
140 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->na);
141 posned = ftell(inm);
142 continue;
143 }
144 else if(!strncmp(chaine,"RequiredEdges",strlen("RequiredEdges"))) {
145 MMG_FSCANF(inm,"%" MMG5_PRId "",&nreqed);
146 posreqed = ftell(inm);
147 continue;
148 }
149 else if(!strncmp(chaine,"RequiredTriangles",strlen("RequiredTriangles"))) {
150 MMG_FSCANF(inm,"%" MMG5_PRId "",&ntreq);
151 posntreq = ftell(inm);
152 continue;
153 }
154 }
155 }
156 else {
157 bdim = 0;
158 MMG_FREAD(&mesh->ver,MMG5_SW,1,inm);
159 iswp=0;
160 if(mesh->ver==16777216)
161 iswp=1;
162 else if(mesh->ver!=1) {
163 fprintf(stdout,"BAD FILE ENCODING\n");
164 }
165 MMG_FREAD(&mesh->ver,MMG5_SW,1,inm);
166 if(iswp) mesh->ver = MMG5_swapbin(mesh->ver);
167 while(fread(&binch,MMG5_SW,1,inm)!=0 && binch!=54 ) {
168 if(iswp) binch=MMG5_swapbin(binch);
169 if(binch==54) break;
170 if(!bdim && binch==3) { //Dimension
171 MMG_FREAD(&bdim,MMG5_SW,1,inm); //NulPos=>20
172 if(iswp) bdim=MMG5_swapbin(bdim);
173 MMG_FREAD(&bdim,MMG5_SW,1,inm);
174 if(iswp) bdim=MMG5_swapbin(bdim);
175 mesh->dim = bdim;
176 if(bdim!=2) {
177 fprintf(stdout,"BAD MESH DIMENSION : %d\n",mesh->dim);
178 return -1;
179 }
180 continue;
181 } else if(!mesh->np && binch==4) { //Vertices
182 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
183 if(iswp) bpos=MMG5_swapbin(bpos);
184 MMG_FREAD(&mesh->np,MMG5_SW,1,inm);
185 if(iswp) mesh->np=MMG5_SWAPBIN(mesh->np);
186 posnp = ftell(inm);
187 rewind(inm);
188 fseek(inm,bpos,SEEK_SET);
189 continue;
190 } else if(!mesh->nt && binch==6) {//MMG5_Triangles
191 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
192 if(iswp) bpos=MMG5_swapbin(bpos);
193 MMG_FREAD(&mesh->nt,MMG5_SW,1,inm);
194 if(iswp) mesh->nt=MMG5_SWAPBIN(mesh->nt);
195 posnt = ftell(inm);
196 rewind(inm);
197 fseek(inm,bpos,SEEK_SET);
198 continue;
199 }
200 else if(binch==17) { //RequiredTriangles
201 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
202 if(iswp) bpos=MMG5_swapbin(bpos);
203 MMG_FREAD(&ntreq,MMG5_SW,1,inm);
204 if(iswp) ntreq=MMG5_SWAPBIN(ntreq);
205 posntreq = ftell(inm);
206 rewind(inm);
207 fseek(inm,bpos,SEEK_SET);
208 continue;
209 } else if(!mesh->nquad && binch==7) {//Quadrilaterals
210 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
211 if(iswp) bpos=MMG5_swapbin(bpos);
212 MMG_FREAD(&mesh->nquad,MMG5_SW,1,inm);
213 if(iswp) mesh->nquad=MMG5_SWAPBIN(mesh->nquad);
214 posnq = ftell(inm);
215 rewind(inm);
216 fseek(inm,bpos,SEEK_SET);
217 continue;
218 } else if(binch==18) { //RequiredQuadrilaterals
219 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
220 if(iswp) bpos=MMG5_swapbin(bpos);
221 MMG_FREAD(&nqreq,MMG5_SW,1,inm);
222 if(iswp) nqreq=MMG5_SWAPBIN(nqreq);
223 posnqreq = ftell(inm);
224 rewind(inm);
225 fseek(inm,bpos,SEEK_SET);
226 continue;
227 } else if(!ncor && binch==13) {
228 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
229 if(iswp) bpos=MMG5_swapbin(bpos);
230 MMG_FREAD(&ncor,MMG5_SW,1,inm);
231 if(iswp) ncor=MMG5_SWAPBIN(ncor);
232 posncor = ftell(inm);
233 rewind(inm);
234 fseek(inm,bpos,SEEK_SET);
235 continue;
236 } else if(!mesh->na && binch==5) { //Edges
237 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
238 if(iswp) bpos=MMG5_swapbin(bpos);
239 MMG_FREAD(&mesh->na,MMG5_SW,1,inm);
240 if(iswp) mesh->na=MMG5_SWAPBIN(mesh->na);
241 posned = ftell(inm);
242 rewind(inm);
243 fseek(inm,bpos,SEEK_SET);
244 continue;
245 } else if(!nreqed && binch==16) { //RequiredEdges
246 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
247 if(iswp) bpos=MMG5_swapbin(bpos);
248 MMG_FREAD(&nreqed,MMG5_SW,1,inm);
249 if(iswp) nreqed=MMG5_SWAPBIN(nreqed);
250 posreqed = ftell(inm);
251 rewind(inm);
252 fseek(inm,bpos,SEEK_SET);
253 continue;
254 } else if(!nreq && binch==15) { //RequiredVertices
255 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
256 if(iswp) bpos=MMG5_swapbin(bpos);
257 MMG_FREAD(&nreq,MMG5_SW,1,inm);
258 if(iswp) nreq=MMG5_SWAPBIN(nreq);
259 posreq = ftell(inm);
260 rewind(inm);
261 fseek(inm,bpos,SEEK_SET);
262 continue;
263 } else {
264 //printf("on traite ? %d\n",binch);
265 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
266 if(iswp) bpos=MMG5_swapbin(bpos);
267 //printf("on avance... Nulpos %d\n",bpos);
268 rewind(inm);
269 fseek(inm,bpos,SEEK_SET);
270 }
271 }
272
273 }
274
275 if ( !mesh->np ) {
276 fprintf(stdout," ** MISSING DATA : no point\n");
277 return -1;
278 }
279
280 mesh->npi = mesh->np;
281 mesh->nai = mesh->na;
282 mesh->nti = mesh->nt;
283 if ( !mesh->np ) {
284 fprintf(stdout," ** MISSING DATA\n");
285 return -1;
286 }
287
288 /* Memory allocation */
289 if ( !MMG2D_zaldy(mesh) ) return -1;
290
291 /* Read vertices */
292 rewind(inm);
293 fseek(inm,posnp,SEEK_SET);
294 for (k=1; k<=mesh->np; k++) {
295 ppt = &mesh->point[k];
296 if (mesh->ver < 2) { /*float*/
297 if (!bin) {
298 if ( mesh->info.renum >=2 ) {
299 for (i=0 ; i<3 ; i++) {
300 MMG_FSCANF(inm,"%f",&fc);
301 if(i==2) break;
302 ppt->c[i] = (double) fc;
303 }
304 } else {
305 for (i=0 ; i<2 ; i++) {
306 MMG_FSCANF(inm,"%f",&fc);
307 ppt->c[i] = (double) fc;
308 }
309 }
310 MMG_FSCANF(inm,"%" MMG5_PRId "",&ppt->ref);
311 } else {
312 if ( mesh->info.renum >= 2 ) {
313 fprintf(stderr," ## Warning: %s: binary not available with"
314 " -msh option.\n",__func__);
315 return -1;
316 }
317 for (i=0 ; i<2 ; i++) {
318 MMG_FREAD(&fc,MMG5_SW,1,inm);
319 if(iswp) fc=MMG5_swapf(fc);
320 ppt->c[i] = (double) fc;
321 }
322 MMG_FREAD(&ppt->ref,MMG5_SW,1,inm);
323 if(iswp) ppt->ref=MMG5_SWAPBIN(ppt->ref);
324 }
325 } else {
326 if (!bin) {
327 if ( mesh->info.renum >= 2 ) {
328 MMG_FSCANF(inm,"%lf %lf %lf %" MMG5_PRId "",&ppt->c[0],&ppt->c[1],&dtmp,&ppt->ref);
329 } else {
330 MMG_FSCANF(inm,"%lf %lf %" MMG5_PRId "",&ppt->c[0],&ppt->c[1],&ppt->ref);
331 }
332 }
333 else {
334 for (i=0 ; i<2 ; i++) {
335 MMG_FREAD(&ppt->c[i],MMG5_SD,1,inm);
336 if(iswp) ppt->c[i]=MMG5_swapd(ppt->c[i]);
337 }
338 MMG_FREAD(&ppt->ref,MMG5_SW,1,inm);
339 if(iswp) ppt->ref=MMG5_SWAPBIN(ppt->ref);
340 }
341 }
342 if ( ppt->ref < 0 ) {
343 ppt->ref = -ppt->ref;
344 ++nref;
345 }
346 ppt->tag = MG_NUL;
347 }
348
349 /* Read edges */
350 rewind(inm);
351 fseek(inm,posned,SEEK_SET);
352 for (k=1; k<=mesh->na; k++) {
353 ped = &mesh->edge[k];
354 if (!bin) {
355 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&ped->a,&ped->b,&ped->ref);
356 }
357 else {
358 MMG_FREAD(&ped->a,MMG5_SW,1,inm);
359 if(iswp) ped->a=MMG5_SWAPBIN(ped->a);
360 MMG_FREAD(&ped->b,MMG5_SW,1,inm);
361 if(iswp) ped->b=MMG5_SWAPBIN(ped->b);
362 MMG_FREAD(&ped->ref,MMG5_SW,1,inm);
363 if(iswp) ped->ref=MMG5_SWAPBIN(ped->ref);
364 }
365 if ( ped->ref < 0 ) {
366 ped->ref = -ped->ref;
367 ++nref;
368 }
369 ped->tag |= MG_REF+MG_BDY;
370 }
371
372 /* Read triangles */
373 if ( mesh->nt ) {
374 rewind(inm);
375 fseek(inm,posnt,SEEK_SET);
376 norient = 0;
377 for (k=1; k<=mesh->nt; k++) {
378 pt = &mesh->tria[k];
379 if (!bin) {
380 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&pt->v[0],&pt->v[1],&pt->v[2],&pt->ref);
381 }
382 else {
383 for (i=0 ; i<3 ; i++) {
384 MMG_FREAD(&pt->v[i],MMG5_SW,1,inm);
385 if(iswp) pt->v[i]=MMG5_SWAPBIN(pt->v[i]);
386 }
387 MMG_FREAD(&pt->ref,MMG5_SW,1,inm);
388 if(iswp) pt->ref=MMG5_SWAPBIN(pt->ref);
389 }
390 for (i=0; i<3; i++) {
391 ppt = &mesh->point[ pt->v[i] ];
392 ppt->tag &= ~MG_NUL;
393 }
394 for(i=0 ; i<3 ; i++)
395 pt->edg[i] = 0;
396
397 /* Get positive ref */
398 if ( pt->ref < 0 ) {
399 pt->ref = -pt->ref;
400 ++nref;
401 }
402
403 /* Check orientation */
404 air = MMG2D_quickarea(mesh->point[pt->v[0]].c,mesh->point[pt->v[1]].c,
405 mesh->point[pt->v[2]].c);
406 if(air < 0) {
407 norient++;
408 tmp = pt->v[2];
409 pt->v[2] = pt->v[1];
410 pt->v[1] = tmp;
411 }
412 }
413 if ( norient ) {
414 fprintf(stdout,"\n $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n");
415 fprintf(stdout," BAD ORIENTATION : vol < 0 -- %8" MMG5_PRId " element(s) reoriented\n",norient);
416 fprintf(stdout," $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n\n");
417 }
418
419 if ( ntreq ) {
420 rewind(inm);
421 fseek(inm,posntreq,SEEK_SET);
422 for (k=1; k<=ntreq; k++) {
423 if(!bin) {
424 MMG_FSCANF(inm,"%" MMG5_PRId "",&i);
425 }
426 else {
427 MMG_FREAD(&i,MMG5_SW,1,inm);
428 if(iswp) i=MMG5_SWAPBIN(i);
429 }
430 if ( i>mesh->nt ) {
431 fprintf(stderr,"\n ## Warning: %s: required triangle number %8" MMG5_PRId ""
432 " ignored.\n",__func__,i);
433 } else {
434 pt = &mesh->tria[i];
435 pt->tag[0] |= MG_REQ;
436 pt->tag[1] |= MG_REQ;
437 pt->tag[2] |= MG_REQ;
438 }
439 }
440 }
441 }
442 else {
443 for (k=1; k<=mesh->np; k++) {
444 ppt = &mesh->point[ k ];
445 ppt->tag &= ~MG_NUL;
446 }
447 }
448
449 /* read mesh quadrilaterals */
450 if ( mesh->nquad ) {
451 rewind(inm);
452 fseek(inm,posnq,SEEK_SET);
453
454 for (k=1; k<=mesh->nquad; k++) {
455 pq1 = &mesh->quadra[k];
456 if (!bin) {
457 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&pq1->v[0],&pq1->v[1],&pq1->v[2],
458 &pq1->v[3],&pq1->ref);
459 }
460 else {
461 for (i=0 ; i<4 ; i++) {
462 MMG_FREAD(&pq1->v[i],MMG5_SW,1,inm);
463 if(iswp) pq1->v[i]=MMG5_SWAPBIN(pq1->v[i]);
464 }
465 MMG_FREAD(&pq1->ref,MMG5_SW,1,inm);
466 if(iswp) pq1->ref=MMG5_SWAPBIN(pq1->ref);
467 }
468 for (i=0; i<4; i++) {
469 ppt = &mesh->point[ pq1->v[i] ];
470 ppt->tag &= ~MG_NUL;
471 }
472
473 if ( pq1->ref < 0 ) {
474 pq1->ref = -pq1->ref;
475 ++nref;
476 }
477 }
478 /* get required quadrilaterals */
479 if(nqreq) {
480 rewind(inm);
481 fseek(inm,posnqreq,SEEK_SET);
482 for (k=1; k<=nqreq; k++) {
483 if(!bin) {
484 MMG_FSCANF(inm,"%" MMG5_PRId "",&i);
485 }
486 else {
487 MMG_FREAD(&i,MMG5_SW,1,inm);
488 if(iswp) i=MMG5_SWAPBIN(i);
489 }
490 if ( i>mesh->nquad ) {
491 fprintf(stderr,"\n ## Warning: %s: required quadrilaterals number"
492 " %8" MMG5_PRId " ignored.\n",__func__,i);
493 } else {
494 pq1 = &mesh->quadra[i];
495 pq1->tag[0] |= MG_REQ;
496 pq1->tag[1] |= MG_REQ;
497 pq1->tag[2] |= MG_REQ;
498 pq1->tag[3] |= MG_REQ;
499 }
500 }
501 }
502 }
503
504 /* Read corners */
505 if ( ncor ) {
506 rewind(inm);
507 fseek(inm,posncor,SEEK_SET);
508 for (k=1; k<=ncor; k++) {
509 if (!bin) {
510 MMG_FSCANF(inm,"%" MMG5_PRId "",&ref);
511 }
512 else {
513 MMG_FREAD(&ref,MMG5_SW,1,inm);
514 if(iswp) ref=MMG5_SWAPBIN(ref);
515 }
516 ppt = &mesh->point[ref];
517 ppt->tag |= MG_CRN;
518 }
519 }
520
521 /* Read required vertices*/
522 if (nreq) {
523 rewind(inm);
524 fseek(inm,posreq,SEEK_SET);
525 for (k=1; k<=nreq; k++) {
526 if (!bin) {
527 MMG_FSCANF(inm,"%" MMG5_PRId "",&ref);
528 }
529 else {
530 MMG_FREAD(&ref,MMG5_SW,1,inm);
531 if(iswp) ref=MMG5_SWAPBIN(ref);
532 }
533 ppt = &mesh->point[ref];
534 ppt->tag |= MG_REQ;
535 ppt->tag &= ~MG_NUL;
536 }
537 }
538
539 /* read required edges*/
540 if (nreqed) {
541 rewind(inm);
542 fseek(inm,posreqed,SEEK_SET);
543 for (k=1; k<=nreqed; k++) {
544 if (!bin) {
545 MMG_FSCANF(inm,"%" MMG5_PRId "",&ref);
546 }
547 else {
548 MMG_FREAD(&ref,MMG5_SW,1,inm);
549 if(iswp) ref=MMG5_SWAPBIN(ref);
550 }
551 ped = &mesh->edge[ref];
552 ped->tag |= MG_REQ;
553 }
554 }
555
556 fclose(inm);
557
558 if ( nref ) {
559 fprintf(stdout,"\n $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n");
560 fprintf(stdout," WARNING : %" MMG5_PRId " entities with unexpected refs (ref< 0).\n",nref);
561 fprintf(stdout," We take their absolute values.\n");
562 fprintf(stdout," $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n\n");
563 }
564
565 if ( abs(mesh->info.imprim) > 3 ) {
566 fprintf(stdout," NUMBER OF VERTICES %8" MMG5_PRId " CORNERS %6" MMG5_PRId "\n",mesh->np,ncor);
567 fprintf(stdout," NUMBER OF TRIANGLES %8" MMG5_PRId "\n",mesh->nt);
568 if ( mesh->nquad )
569 fprintf(stdout," NUMBER OF QUADRILATERALS %8" MMG5_PRId "\n",mesh->nquad);
570
571 if ( mesh->na )
572 fprintf(stdout," NUMBER OF EDGES %8" MMG5_PRId "\n",mesh->na);
573
574 if ( nreq || nreqed || ntreq || nqreq ) {
575 fprintf(stdout," NUMBER OF REQUIRED ENTITIES: \n");
576 if ( nreq )
577 fprintf(stdout," VERTICES %8" MMG5_PRId " \n",nreq);
578 if ( nreqed )
579 fprintf(stdout," EDGES %8" MMG5_PRId " \n",nreqed);
580 if ( ntreq )
581 fprintf(stdout," TRIANGLES %8" MMG5_PRId " \n",ntreq);
582 if ( nqreq )
583 fprintf(stdout," QUADRILATERALS %8" MMG5_PRId " \n",nqreq);
584 }
585 if(ncor)
586 fprintf(stdout," NUMBER OF CORNERS %8" MMG5_PRId " \n",ncor);
587 }
588
589 return 1;
590}
591
593 int ier=0;
594 const char *filenameptr,*solnameptr;
595 char *tmp,*soltmp;
596
597 if ( filename && strlen(filename) ) {
598 filenameptr = filename;
599 solnameptr = filename;
600 }
601 else if (mesh->namein && strlen(mesh->namein) ) {
602 filenameptr = mesh->namein;
603 if ( sol && strlen(sol->namein) ) {
604 solnameptr = sol->namein;
605 }
606 else {
607 solnameptr = mesh->namein;
608 }
609 }
610 else {
611 fprintf(stderr," ## Error: %s: please provide input file name"
612 " (either in the mesh structure or as function argument).\n",
613 __func__);
614 return -1;
615 }
616
617 MMG5_SAFE_MALLOC(tmp,strlen(filenameptr)+1,char,return -1);
618 strcpy(tmp,filenameptr);
619
620 /* read mesh/sol files */
621 char *ptr = MMG5_Get_filenameExt(tmp);
622 int fmtin = MMG5_Get_format(ptr,MMG5_FMT_MeditASCII);
623
624 switch ( fmtin ) {
625
626 case ( MMG5_FMT_GmshASCII ): case ( MMG5_FMT_GmshBinary ):
628 break;
629
630 case ( MMG5_FMT_VtkVtu ):
632 break;
633
634 case ( MMG5_FMT_VtkVtk ):
636 break;
637
638 case ( MMG5_FMT_MeditASCII ): case ( MMG5_FMT_MeditBinary ):
640 if ( ier < 1 ) { break; }
641
642 /* Facultative metric */
643 if ( sol ) {
644 MMG5_SAFE_MALLOC(soltmp,strlen(solnameptr)+1,char,return -1);
645 strcpy(soltmp,solnameptr);
646
647 if ( MMG2D_loadSol(mesh,sol,tmp) == -1) {
648 fprintf(stderr,"\n ## ERROR: WRONG DATA TYPE OR WRONG SOLUTION NUMBER.\n");
649 ier = 0;
650 }
651 MMG5_SAFE_FREE(soltmp);
652 }
653 break;
654
655 default:
656 fprintf(stderr," ** I/O AT FORMAT %s NOT IMPLEMENTED.\n",MMG5_Get_formatName(fmtin) );
657 ier= 0;
658 }
659
661
662 return ier;
663}
664
665
675 MMG5_pPoint ppt;
676 double z;
677 MMG5_int k;
678
679 if (!mesh->nt) {
680 for (k=1; k<=mesh->np; k++) {
681 ppt = &mesh->point[ k ];
682 ppt->tag &= ~MG_NUL;
683 }
684 }
685
686 z = 0.;
687 for ( k=1; k<=mesh->np; ++k ) {
688 ppt = &mesh->point[k];
689 if ( !MG_VOK(ppt) ) continue;
690
691 z += fabs(ppt->c[2]);
692 }
693 if ( z > MMG5_EPSOK ) {
694 fprintf(stderr,"\n ## Error: %s: Input mesh must be a two-dimensional mesh.\n",
695 __func__);
696 return 0;
697 }
698 return 1;
699}
700
702 FILE* inm;
703 long posNodes,posElts,*posNodeData;
704 int ier;
705 int bin,iswp,nsols;
706 MMG5_int nelts;
707
708 mesh->dim = 2;
709
711 &posNodes,&posElts,&posNodeData,
712 &bin,&iswp,&nelts,&nsols);
713 if ( ier < 1 ) return (ier);
714
715 if ( nsols>1 ) {
716 fprintf(stderr,"Error: SEVERAL SOLUTIONS FOUND (%d)\n",nsols);
717 fclose(inm);
718 MMG5_SAFE_FREE(posNodeData);
719 return -1;
720 }
721
722 if ( !MMG2D_zaldy(mesh) ) {
723 fclose(inm);
724 MMG5_SAFE_FREE(posNodeData);
725 return -1;
726 }
727
728 if ( mesh->ne || mesh->nprism ) {
729 fprintf(stderr,"\n ## Error: %s: Input mesh must be a two-dimensional mesh.\n",
730 __func__);
731 fclose(inm);
732 MMG5_SAFE_FREE(posNodeData);
733 return -1;
734 }
735 if ( !mesh->nt )
736 fprintf(stdout," ** WARNING NO GIVEN TRIANGLE\n");
737
738 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt ) {
739 fclose(inm);
740 MMG5_SAFE_FREE(posNodeData);
741 return -1;
742 }
743
745 posNodes,posElts,posNodeData,
746 bin,iswp,nelts,nsols);
747
748 MMG5_SAFE_FREE(posNodeData);
749 if ( ier < 1 ) return ier;
750
751 if ( sol ) {
752 /* Check the metric type */
753 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,inm);
754 if ( ier <1 ) {
755 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
756 return ier;
757 }
758 }
759
760 /* Mark all points as used in case of mesh generation and check the
761 * z-componant */
762 if ( !MMG2D_2dMeshCheck(mesh) ) return -1;
763
764 return 1;
765}
766
768 FILE* inm;
769 long posNodes,posElts,*posNodeData;
770 int ier;
771 int bin,iswp,nsols;
772 MMG5_int nelts;
773
774 mesh->dim = 2;
775
777 &posNodes,&posElts,&posNodeData,
778 &bin,&iswp,&nelts,&nsols);
779 if ( ier < 1 ) return (ier);
780
781 mesh->nsols = nsols;
782 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
783
784 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
785 printf(" Exit program.\n"); fclose(inm);
786 MMG5_SAFE_FREE(posNodeData);
787 return -1);
788 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
789
790 if ( !MMG2D_zaldy(mesh) ) {
791 fclose(inm);
792 MMG5_SAFE_FREE(posNodeData);
793 return -1;
794 }
795
796 if ( mesh->ne || mesh->nprism ) {
797 fprintf(stderr,"\n ## Error: %s: Input mesh must be a two-dimensional mesh.\n",
798 __func__);
799 fclose(inm);
800 MMG5_SAFE_FREE(posNodeData);
801 return -1;
802 }
803 if ( !mesh->nt )
804 fprintf(stdout," ** WARNING NO GIVEN TRIANGLE\n");
805
806 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt ) {
807 fclose(inm);
808 MMG5_SAFE_FREE(posNodeData);
809 return -1;
810 }
811
813 posNodes,posElts,posNodeData,
814 bin,iswp,nelts,nsols);
815
816 MMG5_SAFE_FREE(posNodeData);
817 if ( ier < 1 ) {
818 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
819 return ier;
820 }
821
822 /* Mark all points as used in case of mesh generation and check the
823 * z-componant */
824 if ( !MMG2D_2dMeshCheck(mesh) ) return -1;
825
826 return 1;
827}
828
841static inline
842int MMG2D_readFloatSol(MMG5_pSol sol,FILE *inm,int bin,int iswp,MMG5_int pos) {
843 float fbuf;
844 int i;
845
846 for (i=0; i<sol->size; i++) {
847 if ( !bin ) {
848 MMG_FSCANF(inm,"%f",&fbuf);
849 sol->m[sol->size*pos+i] = (double)fbuf;
850 }
851 else {
852 MMG_FREAD(&fbuf,MMG5_SW,1,inm);
853 if ( iswp ) fbuf=MMG5_swapf(fbuf);
854 sol->m[sol->size*pos+i] = (double)fbuf;
855 }
856 }
857 return 1;
858}
859
872static inline
873int MMG2D_readDoubleSol(MMG5_pSol sol,FILE *inm,int bin,int iswp,MMG5_int pos) {
874 double dbuf;
875 int i;
876
877 for (i=0; i<sol->size; i++) {
878 if ( !bin ) {
879 MMG_FSCANF(inm,"%lf",&dbuf);
880 sol->m[sol->size*pos+i] = (double)dbuf;
881 }
882 else {
883 MMG_FREAD(&dbuf,MMG5_SD,1,inm);
884 if ( iswp ) dbuf=MMG5_swapf(dbuf);
885 sol->m[sol->size*pos+i] = (double)dbuf;
886 }
887 }
888 return 1;
889}
890
901 FILE *inm;
902 long posnp;
903 int iswp,ier,meshDim,*type,nsols,dim;
904 int ver,bin;
905 MMG5_int k,np;
906
908 meshDim = 2;
909 if ( mesh->info.renum >= 2 ) {
910 /* -msh mode */
911 meshDim = 3;
912 }
913 ier = MMG5_loadSolHeader(filename,meshDim,&inm,&ver,&bin,&iswp,&np,&dim,&nsols,
914 &type,&posnp,mesh->info.imprim);
915
916 /* correction for the -msh mode */
917 sol->dim = 2;
918
919 if ( ier < 1 ) return ier;
920
921 if ( nsols!=1 ) {
922 fprintf(stderr,"Error: SEVERAL SOLUTIONS FOUND (%d)\n",nsols);
923 fclose(inm);
924 MMG5_SAFE_FREE(type);
925 return -1;
926 }
927
928 if ( mesh->np != np ) {
929 fprintf(stderr," ** MISMATCHES DATA: THE NUMBER OF VERTICES IN "
930 "THE MESH (%" MMG5_PRId ") DIFFERS FROM THE NUMBER OF VERTICES IN "
931 "THE SOLUTION (%" MMG5_PRId ") \n",mesh->np,np);
932 fclose(inm);
933 MMG5_SAFE_FREE(type);
934 return -1;
935 }
936
937 /* #MMG5_loadSolHeader function reads only solution at vertices so we don't
938 have to check the entites on which the metric applies */
939 int entities = MMG5_Vertex;
940 ier = MMG5_chkMetricType(mesh,type,&entities,inm);
941 if ( ier < 1 ) {
942 MMG5_SAFE_FREE(type);
943 return ier;
944 }
945
946 /* Allocate and store the header informations for each solution */
947 if ( !MMG2D_Set_solSize(mesh,sol,MMG5_Vertex,mesh->np,type[0]) ) {
948 fclose(inm);
949 MMG5_SAFE_FREE(type);
950 return -1;
951 }
952 /* For binary file, we read the verson inside the file */
953 if ( ver ) sol->ver = ver;
954
955 MMG5_SAFE_FREE(type);
956
957 /* Read mesh solutions */
958 rewind(inm);
959 fseek(inm,posnp,SEEK_SET);
960
961 if ( sol->ver == 1 ) {
962 /* Simple precision */
963 for (k=1; k<=sol->np; k++) {
964 if ( MMG2D_readFloatSol(sol,inm,bin,iswp,k) < 0 ) return -1;
965 }
966 }
967 else {
968 for (k=1; k<=sol->np; k++) {
969 /* Double precision */
970 if ( MMG2D_readDoubleSol(sol,inm,bin,iswp,k) < 0 ) return -1;
971 }
972 }
973
974 fclose(inm);
975
976 /* stats */
978
979 return 1;
980}
981
992 MMG5_pSol psl;
993 FILE *inm;
994 long posnp;
995 int iswp,ier,meshDim,nsols,*type;
996 MMG5_int k,np;
997 int j,ver,bin,dim;
998 char data[16];
999 static int8_t mmgWarn = 0;
1000
1002 meshDim = 2;
1003 if ( mesh->info.renum >= 2 ) {
1004 /* -msh mode */
1005 meshDim = 3;
1006 }
1007
1008 ier = MMG5_loadSolHeader(filename,meshDim,&inm,&ver,&bin,&iswp,&np,&dim,&nsols,
1009 &type,&posnp,mesh->info.imprim);
1010 if ( ier < 1 ) return ier;
1011
1012 if ( mesh->np != np ) {
1013 fprintf(stderr," ** MISMATCHES DATA: THE NUMBER OF VERTICES IN "
1014 "THE MESH (%" MMG5_PRId ") DIFFERS FROM THE NUMBER OF VERTICES IN "
1015 "THE SOLUTION (%" MMG5_PRId ") \n",mesh->np,np);
1016 fclose(inm);
1017 MMG5_SAFE_FREE(type);
1018 return -1;
1019 }
1020
1022 mesh->nsols = nsols;
1023
1024 if ( nsols > MMG5_NSOLS_MAX ) {
1025 fprintf(stderr,"\n ## Error: %s: unexpected number of data (%d).\n",
1026 __func__,nsols);
1027 MMG5_SAFE_FREE(type);
1028 fclose(inm);
1029 return -1;
1030 }
1031
1032 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
1033
1034 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
1035 printf(" Exit program.\n"); fclose(inm);
1036 MMG5_SAFE_FREE(type);
1037 return -1);
1038 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
1039
1040 for ( j=0; j<nsols; ++j) {
1041 psl = *sol+j;
1042
1043 /* Give an arbitrary name to the solution because the Medit format has non
1044 * name field */
1045 sprintf(data,"sol_%d",j);
1046 if ( !MMG2D_Set_inputSolName(mesh,psl,data) ) {
1047 if ( !mmgWarn ) {
1048 mmgWarn = 1;
1049 fprintf(stderr,"\n ## Warning: %s: unable to set solution name for"
1050 " at least 1 solution.\n",__func__);
1051 }
1052 }
1053
1054 /* Allocate and store the header informations for each solution */
1055 if ( !MMG2D_Set_solSize(mesh,psl,MMG5_Vertex,mesh->np,type[j]) ) {
1056 MMG5_SAFE_FREE(type);
1057 fclose(inm);
1058 return -1;
1059 }
1060 psl->dim = 2;
1061 /* For binary file, we read the verson inside the file */
1062 if ( ver ) psl->ver = ver;
1063 }
1064 MMG5_SAFE_FREE(type);
1065
1066 /* read mesh solutions */
1067 rewind(inm);
1068 fseek(inm,posnp,SEEK_SET);
1069
1070 if ( (*sol)[0].ver == 1 ) {
1071 /* Simple precision */
1072 for (k=1; k<=mesh->np; k++) {
1073 for ( j=0; j<nsols; ++j ) {
1074 psl = *sol+j;
1075 if ( MMG2D_readFloatSol(psl,inm,bin,iswp,k) < 0 ) return -1;
1076 }
1077 }
1078 }
1079 else {
1080 /* Double precision */
1081 for (k=1; k<=mesh->np; k++) {
1082 for ( j=0; j<nsols; ++j ) {
1083 psl = *sol+j;
1084 if ( MMG2D_readDoubleSol(psl,inm,bin,iswp,k) < 0 ) return -1;
1085 }
1086 }
1087 }
1088 fclose(inm);
1089
1090 /* stats */
1092
1093 return 1;
1094}
1095
1097 FILE* inm;
1098 MMG5_pPoint ppt;
1099 MMG5_pEdge ped;
1100 MMG5_pTria pt;
1101 MMG5_pQuad pq;
1102 double dblb;
1103 MMG5_int k,ne,np,nc,nreq,nereq,nedreq,nq,nqreq,bpos,ref;
1104 int bin, binch,gmsh;
1105 char *ptr,*data,chaine[MMG5_FILESTR_LGTH];
1106
1107 gmsh = (mesh->info.renum==1||mesh->info.renum==2);
1108
1109 mesh->ver = 2;
1110 bin = 0;
1111
1112 /* Name of file */
1113 MMG5_SAFE_CALLOC(data,strlen(filename)+7,char,return 0);
1114 strcpy(data,filename);
1115 ptr = strstr(data,".mesh");
1116 if ( !ptr ) {
1117 strcat(data,".meshb");
1118 if( !(inm = fopen(data,"wb")) ) {
1119 ptr = strstr(data,".mesh");
1120 *ptr = '\0';
1121 strcat(data,".mesh");
1122 if( !(inm = fopen(data,"wb")) ) {
1123 MMG5_SAFE_FREE(data);
1124 return 0;
1125 }
1126 }
1127 else {
1128 bin = 1;
1129 }
1130 }
1131 else {
1132 ptr = strstr(data,".meshb");
1133 if( ptr ) bin = 1;
1134 if( !(inm = fopen(data,"wb")) ) {
1135 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
1136 MMG5_SAFE_FREE(data);
1137 return 0;
1138 }
1139 }
1140
1141 if ( mesh->info.imprim >= 0 )
1142 fprintf(stdout," %%%% %s OPENED\n",data);
1143 MMG5_SAFE_FREE(data);
1144
1145 /* Write header */
1146 binch=0; bpos=10;
1147 if ( !bin ) {
1148 strcpy(&chaine[0],"MeshVersionFormatted 2\n");
1149 fprintf(inm,"%s",chaine);
1150 if ( gmsh ) {
1151 strcpy(&chaine[0],"\n\nDimension 3\n");
1152 }
1153 else {
1154 strcpy(&chaine[0],"\n\nDimension 2\n");
1155 }
1156 fprintf(inm,"%s ",chaine);
1157 }
1158 else
1159 {
1160 binch = 1; //MeshVersionFormatted
1161 fwrite(&binch,MMG5_SW,1,inm);
1162 binch = 2; //version
1163 fwrite(&binch,MMG5_SW,1,inm);
1164 binch = 3; //Dimension
1165 fwrite(&binch,MMG5_SW,1,inm);
1166 bpos = 20; //Pos
1167 fwrite(&bpos,MMG5_SW,1,inm);
1168 if ( gmsh ) {
1169 binch = 3; //Dimension
1170 }
1171 else binch = 2;
1172 fwrite(&binch,MMG5_SW,1,inm);
1173 }
1174
1175 /* Write vertices */
1176 np = 0;
1177 for (k=1; k<=mesh->np; k++) {
1178 ppt = &mesh->point[k];
1179 if ( MG_VOK(ppt) ) np++;
1180 ppt->tmp = np;
1181 }
1182
1183 if ( !bin ) {
1184 strcpy(&chaine[0],"\n\nVertices\n");
1185 fprintf(inm,"%s",chaine);
1186 fprintf(inm,"%" MMG5_PRId "\n",np);
1187 }
1188 else {
1189 binch = 4; //Vertices
1190 fwrite(&binch,MMG5_SW,1,inm);
1191 if ( gmsh )
1192 bpos += (3+(1+3*mesh->ver)*np)*MMG5_SW; //NullPos
1193 else
1194 bpos += (3+(1+2*mesh->ver)*np)*MMG5_SW; //NullPos
1195
1196 fwrite(&bpos,MMG5_SW,1,inm);
1197 fwrite(&np,MMG5_SW,1,inm);
1198 }
1199 fflush(inm);
1200
1201 for (k=1; k<=mesh->np; k++) {
1202 ppt = &mesh->point[k];
1203 if ( MG_VOK(ppt) ) {
1204 ref = ppt->ref;
1205 if ( gmsh ) {
1206 if ( !bin )
1207 fprintf(inm,"%.15lg %.15lg 0 %" MMG5_PRId "\n",ppt->c[0],ppt->c[1],ref);
1208 else {
1209 dblb = 0.;
1210 fwrite((unsigned char*)&ppt->c[0],MMG5_SD,1,inm);
1211 fwrite((unsigned char*)&ppt->c[1],MMG5_SD,1,inm);
1212 fwrite((unsigned char*)&dblb,MMG5_SD,1,inm);
1213 fwrite((unsigned char*)&ref,MMG5_SW,1,inm);
1214 }
1215 }
1216 else {
1217 if ( !bin ) {
1218 fprintf(inm,"%.15lg %.15lg %" MMG5_PRId "\n",ppt->c[0],ppt->c[1],ref);
1219 fflush(inm);
1220 }
1221 else {
1222 fwrite(&ppt->c[0],MMG5_SD,1,inm);
1223 fwrite(&ppt->c[1],MMG5_SD,1,inm);
1224 fwrite(&ref,MMG5_SW,1,inm);
1225 }
1226 }
1227 }
1228 }
1229
1230 /* Print corners */
1231 nc = 0;
1232 for (k=1; k<=mesh->np; k++) {
1233 ppt = &mesh->point[k];
1234 if ( MG_VOK(ppt) && (ppt->tag & MG_CRN) ) nc++;
1235 }
1236
1237 if ( nc ) {
1238 if ( !bin ) {
1239 strcpy(&chaine[0],"\n\nCorners\n");
1240 fprintf(inm,"%s",chaine);
1241 fprintf(inm,"%" MMG5_PRId "\n",nc);
1242 }
1243 else
1244 {
1245 binch = 13; //
1246 fwrite(&binch,MMG5_SW,1,inm);
1247 bpos += (3+nc)*MMG5_SW; //NullPos
1248 fwrite(&bpos,MMG5_SW,1,inm);
1249 fwrite(&nc,MMG5_SW,1,inm);
1250 }
1251
1252 for (k=1; k<=mesh->np; k++) {
1253 ppt = &mesh->point[k];
1254 if ( MG_VOK(ppt) && (ppt->tag & MG_CRN) ) {
1255 if(!bin) {
1256 fprintf(inm,"%" MMG5_PRId "\n",ppt->tmp);
1257 }
1258 else {
1259 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1260 }
1261 }
1262 }
1263 }
1264
1265 /* Required vertex */
1266 nreq = 0;
1267 for (k=1; k<=mesh->np; k++) {
1268 ppt = &mesh->point[k];
1269 if ( MG_VOK(ppt) ) {
1270 if ( mesh->info.nosurf && (ppt->tag & MG_NOSURF) ) continue;
1271 if ( ppt->tag & MG_REQ ) nreq++;
1272 }
1273 }
1274 if ( nreq ) {
1275 if ( !bin ) {
1276 strcpy(&chaine[0],"\n\nRequiredVertices\n");
1277 fprintf(inm,"%s",chaine);
1278 fprintf(inm,"%" MMG5_PRId "\n",nreq);
1279 }
1280 else {
1281 binch = 15; //
1282 fwrite(&binch,MMG5_SW,1,inm);
1283 bpos += (3+nreq)*MMG5_SW; //NullPos
1284 fwrite(&bpos,MMG5_SW,1,inm);
1285 fwrite(&nreq,MMG5_SW,1,inm);
1286 }
1287 for (k=1; k<=mesh->np; k++) {
1288 ppt = &mesh->point[k];
1289 if ( MG_VOK(ppt) ) {
1290 if ( mesh->info.nosurf && ( ppt->tag & MG_NOSURF )) continue;
1291 if ((ppt->tag & MG_REQ)
1292 /*&& ( (ppt->tag & MG_BDY) || (ppt->tag & MG_SD) ) */ ) {
1293 if(!bin)
1294 fprintf(inm,"%" MMG5_PRId "\n",ppt->tmp);
1295 else
1296 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1297 }
1298 }
1299 }
1300 }
1301
1302 /* edges */
1303 nedreq = 0;
1304 if ( mesh->na ) {
1305 if(!bin) {
1306 strcpy(&chaine[0],"\n\nEdges\n");
1307 fprintf(inm,"%s",chaine);
1308 fprintf(inm,"%" MMG5_PRId "\n",mesh->na);
1309 }
1310 else {
1311 binch = 5; //Edges
1312 fwrite(&binch,MMG5_SW,1,inm);
1313 bpos += (3+3*mesh->na)*MMG5_SW;//Pos
1314 fwrite(&bpos,MMG5_SW,1,inm);
1315 fwrite(&mesh->na,MMG5_SW,1,inm);
1316 }
1317 for (k=1; k<=mesh->na; k++) {
1318 ped = &mesh->edge[k];
1319 if(!bin)
1320 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",mesh->point[ped->a].tmp,mesh->point[ped->b].tmp,ped->ref);
1321 else
1322 {
1323 fwrite(&mesh->point[ped->a].tmp,MMG5_SW,1,inm);
1324 fwrite(&mesh->point[ped->b].tmp,MMG5_SW,1,inm);
1325 fwrite(&ped->ref,MMG5_SW,1,inm);
1326 }
1327 if ( ped->tag & MG_REQ ) nedreq++;
1328 }
1329
1330 if ( nedreq ) {
1331 if(!bin) {
1332 strcpy(&chaine[0],"\n\nRequiredEdges\n");
1333 fprintf(inm,"%s",chaine);
1334 fprintf(inm,"%" MMG5_PRId "\n",nedreq);
1335 } else {
1336 binch = 16; //RequiredEdges
1337 fwrite(&binch,MMG5_SW,1,inm);
1338 bpos += (3 + nedreq)*MMG5_SW;//Pos
1339 fwrite(&bpos,MMG5_SW,1,inm);
1340 fwrite(&nedreq,MMG5_SW,1,inm);
1341 }
1342 ne = 0;
1343 for (k=1; k<=mesh->na; k++) {
1344 ne++;
1345 if ( mesh->edge[k].tag & MG_REQ ) {
1346 if(!bin) {
1347 fprintf(inm,"%" MMG5_PRId "\n",ne);
1348 } else {
1349 fwrite(&ne,MMG5_SW,1,inm);
1350 }
1351 }
1352 }
1353 }
1354 }
1355
1356 /* elements */
1357 ne = 0;
1358 nereq = 0;
1359 for (k=1; k<=mesh->nt; k++) {
1360 pt = &mesh->tria[k];
1361 if ( !MG_EOK(pt) ) continue;
1362 ne++;
1363 if ( (pt->tag[0] & MG_REQ) && (pt->tag[1] & MG_REQ) && pt->tag[2] & MG_REQ )
1364 ++nereq;
1365 }
1366
1367 if ( ne ) {
1368 if ( !bin ) {
1369 strcpy(&chaine[0],"\n\nTriangles\n");
1370 fprintf(inm,"%s",chaine);
1371 fprintf(inm,"%" MMG5_PRId "\n",ne);
1372 }
1373 else {
1374 binch = 6; //Triangles
1375 fwrite(&binch,MMG5_SW,1,inm);
1376 bpos += (3+4*ne)*MMG5_SW; //Pos
1377 fwrite(&bpos,MMG5_SW,1,inm);
1378 fwrite(&ne,MMG5_SW,1,inm);
1379 }
1380 for (k=1; k<=mesh->nt; k++) {
1381 pt = &mesh->tria[k];
1382 if ( MG_EOK(pt) ) {
1383 ref = pt->ref;
1384 if ( !bin ) {
1385 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",mesh->point[pt->v[0]].tmp,
1386 mesh->point[pt->v[1]].tmp,
1387 mesh->point[pt->v[2]].tmp,ref);
1388 }
1389 else {
1390 fwrite(&mesh->point[pt->v[0]].tmp,MMG5_SW,1,inm);
1391 fwrite(&mesh->point[pt->v[1]].tmp,MMG5_SW,1,inm);
1392 fwrite(&mesh->point[pt->v[2]].tmp,MMG5_SW,1,inm);
1393 fwrite(&ref,MMG5_SW,1,inm);
1394 }
1395 }
1396 }
1397 if ( nereq ) {
1398 if(!bin) {
1399 strcpy(&chaine[0],"\n\nRequiredTriangles\n");
1400 fprintf(inm,"%s",chaine);
1401 fprintf(inm,"%" MMG5_PRId "\n",nereq);
1402 } else {
1403 binch = 17; //ReqTriangles
1404 fwrite(&binch,MMG5_SW,1,inm);
1405 bpos += (3+nereq)*MMG5_SW; //Pos
1406 fwrite(&bpos,MMG5_SW,1,inm);
1407 fwrite(&nereq,MMG5_SW,1,inm);
1408 }
1409 ne = 0;
1410 for (k=1; k<=mesh->nt; k++) {
1411 pt = &mesh->tria[k];
1412 if ( !MG_EOK(pt) ) continue;
1413 ++ne;
1414 if ( (pt->tag[0] & MG_REQ) && (pt->tag[1] & MG_REQ)
1415 && pt->tag[2] & MG_REQ ) {
1416 if(!bin) {
1417 fprintf(inm,"%" MMG5_PRId "\n",ne);
1418 } else {
1419 fwrite(&ne,MMG5_SW,1,inm);
1420 }
1421 }
1422 }
1423 }
1424 }
1425
1426 /* quad + required quad */
1427 nq = nqreq = 0;
1428
1429 if ( mesh->nquad ) {
1430
1431 for (k=1; k<=mesh->nquad; k++) {
1432 pq = &mesh->quadra[k];
1433 if ( !MG_EOK(pq) ) {
1434 continue;
1435 }
1436 nq++;
1437 if ( pq->tag[0] & MG_REQ && pq->tag[1] & MG_REQ &&
1438 pq->tag[2] & MG_REQ && pq->tag[3] & MG_REQ ) {
1439 nqreq++;
1440 }
1441 }
1442 }
1443
1444 if ( nq ) {
1445 if(!bin) {
1446 strcpy(&chaine[0],"\n\nQuadrilaterals\n");
1447 fprintf(inm,"%s",chaine);
1448 fprintf(inm,"%" MMG5_PRId "\n",nq);
1449 } else {
1450 binch = 7; //Quadrilaterals
1451 fwrite(&binch,MMG5_SW,1,inm);
1452 bpos += (3+5*nq)*MMG5_SW; //Pos
1453 fwrite(&bpos,MMG5_SW,1,inm);
1454 fwrite(&nq,MMG5_SW,1,inm);
1455 }
1456 for (k=1; k<=mesh->nquad; k++) {
1457 pq = &mesh->quadra[k];
1458 if ( !MG_EOK(pq) ) continue;
1459
1460 if(!bin) {
1461 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",mesh->point[pq->v[0]].tmp,
1462 mesh->point[pq->v[1]].tmp,mesh->point[pq->v[2]].tmp,
1463 mesh->point[pq->v[3]].tmp, pq->ref);
1464 } else {
1465 fwrite(&mesh->point[pq->v[0]].tmp,MMG5_SW,1,inm);
1466 fwrite(&mesh->point[pq->v[1]].tmp,MMG5_SW,1,inm);
1467 fwrite(&mesh->point[pq->v[2]].tmp,MMG5_SW,1,inm);
1468 fwrite(&mesh->point[pq->v[3]].tmp,MMG5_SW,1,inm);
1469 fwrite(&pq->ref,MMG5_SW,1,inm);
1470 }
1471 if ( pq->tag[0] & MG_REQ && pq->tag[1] & MG_REQ &&
1472 pq->tag[2] & MG_REQ && pq->tag[3] & MG_REQ ) {
1473 nqreq++;
1474 }
1475 }
1476 if ( nqreq ) {
1477 if(!bin) {
1478 strcpy(&chaine[0],"\n\nRequiredQuadrilaterals\n");
1479 fprintf(inm,"%s",chaine);
1480 fprintf(inm,"%" MMG5_PRId "\n",nqreq);
1481 } else {
1482 binch = 18; //ReqQuad
1483 fwrite(&binch,MMG5_SW,1,inm);
1484 bpos += (3+nqreq)*MMG5_SW; //Pos
1485 fwrite(&bpos,MMG5_SW,1,inm);
1486 fwrite(&nqreq,MMG5_SW,1,inm);
1487 }
1488 for (k=0; k<=mesh->nquad; k++) {
1489 pq = &mesh->quadra[k];
1490 if ( (pq->tag[0] & MG_REQ) && (pq->tag[1] & MG_REQ)
1491 && pq->tag[2] & MG_REQ && pq->tag[3] & MG_REQ ) {
1492 if(!bin) {
1493 fprintf(inm,"%" MMG5_PRId "\n",k);
1494 } else {
1495 fwrite(&k,MMG5_SW,1,inm);
1496 }
1497 }
1498 }
1499 }
1500 }
1501
1502 if(!bin) {
1503 strcpy(&chaine[0],"\n\nEnd\n");
1504 fprintf(inm,"%s",chaine);
1505 }
1506 else {
1507 binch = 54; //End
1508 fwrite(&binch,MMG5_SW,1,inm);
1509 bpos += 2*MMG5_SW; //bpos + End key
1510 fwrite(&bpos,MMG5_SW,1,inm);
1511 }
1512
1513 if ( abs(mesh->info.imprim) > 4 ) {
1514 fprintf(stdout," NUMBER OF VERTICES %8" MMG5_PRId " CORNERS %8" MMG5_PRId ""
1515 " REQUIRED %8" MMG5_PRId "\n",np,nc,nreq);
1516
1517 if ( mesh->na )
1518 fprintf(stdout," NUMBER OF EDGES %8" MMG5_PRId " REQUIRED %8" MMG5_PRId "\n",mesh->na,nedreq);
1519 if ( mesh->nt )
1520 fprintf(stdout," NUMBER OF TRIANGLES %8" MMG5_PRId " REQUIRED %8" MMG5_PRId "\n",
1521 mesh->nt, nereq);
1522 if ( nq )
1523 fprintf(stdout," NUMBER OF QUADRILATERALS %8" MMG5_PRId " REQUIRED %8" MMG5_PRId "\n",nq,nqreq);
1524 }
1525
1526 fclose(inm);
1527
1528 return 1;
1529}
1530
1532 return MMG5_saveMshMesh(mesh,&sol,filename,1);
1533}
1534
1537}
1538
1549static inline
1550void MMG2D_writeDoubleSol(MMG5_pSol sol,FILE *inm,int bin,MMG5_int pos,int gmsh) {
1551 int i,isol;
1552
1553 isol = pos * sol->size;
1554
1555 if ( gmsh ) {
1556 if ( !bin ) {
1557 switch ( sol->size ) {
1558 case 1:
1559 fprintf(inm," %.15lg",sol->m[isol]);
1560 break;
1561 case 2:
1562 fprintf(inm," %.15lg %.15lg 0",sol->m[isol],sol->m[isol+1]);
1563 break;
1564 case 3:
1565 fprintf(inm," %.15lg %.15lg %.15lg 0 0 1",sol->m[isol],sol->m[isol+1],sol->m[isol+2]);
1566 break;
1567 }
1568 }
1569 else {
1570 double dbuf = 0.;
1571
1572 switch ( sol->size ) {
1573 case 1:
1574 fwrite(&sol->m[isol],MMG5_SD,1,inm);
1575 break;
1576 case 2:
1577 fwrite(&sol->m[isol],MMG5_SD,2,inm);
1578 fwrite(&dbuf,MMG5_SD,1,inm);
1579 break;
1580 case 3:
1581 fwrite(&sol->m[isol],MMG5_SD,3,inm);
1582 fwrite(&dbuf,MMG5_SD,1,inm);
1583 fwrite(&dbuf,MMG5_SD,1,inm);
1584 dbuf = 1.;
1585 fwrite(&dbuf,MMG5_SD,1,inm);
1586 break;
1587 }
1588 }
1589 }
1590 else {
1591 if ( !bin ) {
1592 for (i=0; i<sol->size; i++)
1593 fprintf(inm," %.15lg",sol->m[isol + i]);
1594 }
1595 else {
1596 for (i=0; i<sol->size; i++)
1597 fwrite(&sol->m[isol + i],MMG5_SD,1,inm);
1598 }
1599 }
1600}
1601
1602
1613 FILE* inm;
1614 MMG5_pPoint ppt;
1615 MMG5_int k,bpos;
1616 int binch,bin,dim,ier,gmsh;
1617
1618 if ( !sol->np ) return 1;
1619
1620
1621 if ( !(sol->np || sol->m) ) {
1622 fprintf(stderr,"\n ## Warning: %s: no metric data to save.\n",__func__);
1623 return 1;
1624 }
1625
1626 gmsh = (mesh->info.renum==1||mesh->info.renum==2);
1627
1628 sol->ver = 2;
1629
1630 if ( sol->dim==2 && gmsh ) {
1631 dim = 3;
1632 }
1633 else {
1634 dim = 2;
1635 }
1636
1637 bpos = 0;
1638 ier = MMG5_saveSolHeader( mesh,filename,&inm,sol->ver,&bin,&bpos,mesh->np,dim,
1639 1,&sol->entities,&sol->type,&sol->size);
1640
1641 if ( ier < 1 ) return ier;
1642
1643 for (k=1; k<=mesh->np; k++) {
1644 ppt = &mesh->point[k];
1645 if ( !MG_VOK(ppt) ) continue;
1646
1647 MMG2D_writeDoubleSol(sol,inm,bin,k,gmsh);
1648 fprintf(inm,"\n");
1649 }
1650
1651 /* End file */
1652 if ( !bin ) {
1653 fprintf(inm,"\n\nEnd\n");
1654 }
1655 else {
1656 binch = 54; //End
1657 fwrite(&binch,MMG5_SW,1,inm);
1658 }
1659 fclose(inm);
1660
1661 return 1;
1662}
1663
1674 FILE* inm;
1675 MMG5_pPoint ppt;
1676 MMG5_pSol psl;
1677 MMG5_int j,k,bpos;
1678 int binch,bin,npointSols,ncellSols,*size;
1679 int *type,*entities,ier,gmsh,dim;
1680
1681
1682 if ( !(*sol)[0].np ) return 1;
1683
1684 MMG5_SAFE_CALLOC(type,mesh->nsols,int,return 0);
1685 MMG5_SAFE_CALLOC(size,mesh->nsols,int,MMG5_SAFE_FREE(type);return 0);
1686 MMG5_SAFE_CALLOC(entities,mesh->nsols,int,
1687 MMG5_SAFE_FREE(type);MMG5_SAFE_FREE(size);return 0);
1688
1689 npointSols = 0;
1690 ncellSols = 0;
1691
1692 gmsh = (mesh->info.renum==1||mesh->info.renum==2);
1693
1694 if ( gmsh ) {
1695 dim = 3;
1696 }
1697 else {
1698 dim = 2;
1699 }
1700
1701 for (k=0; k<mesh->nsols; ++k ) {
1702 (*sol)[k].ver = 2;
1703
1704 if ( ((*sol)[k].entities==MMG5_Noentity) || ((*sol)[k].entities==MMG5_Vertex) ) {
1705 ++npointSols;
1706 }
1707 else if ( (*sol)[k].entities == MMG5_Triangle ) {
1708 ++ncellSols;
1709 }
1710 else {
1711 printf("\n ## Warning: %s: unexpected entity type for solution %" MMG5_PRId ": %s."
1712 "\n Ignored.\n",
1713 __func__,k,MMG5_Get_entitiesName((*sol)[k].entities));
1714 }
1715
1716 type[k] = (*sol)[k].type;
1717 size[k] = (*sol)[k].size;
1718 entities[k] = (*sol)[k].entities;
1719 }
1720
1721 bpos = 0;
1722 ier = MMG5_saveSolHeader( mesh,filename,&inm,(*sol)[0].ver,&bin,&bpos,mesh->np,
1723 dim,mesh->nsols,entities,type,size);
1724
1725
1726 if ( ier < 1 ) return ier;
1727
1728 for (k=1; k<=mesh->np; k++) {
1729 ppt = &mesh->point[k];
1730 if ( !MG_VOK(ppt) ) continue;
1731
1732 for ( j=0; j<mesh->nsols; ++j ) {
1733 psl = *sol + j;
1734
1735 if ( (psl->entities==MMG5_Noentity) || (psl->entities==MMG5_Vertex) ) {
1736 MMG2D_writeDoubleSol(psl,inm,bin,k,gmsh);
1737 }
1738 }
1739 fprintf(inm,"\n");
1740 }
1741
1742 MMG5_saveSolAtTrianglesHeader( mesh,inm,(*sol)[0].ver,bin,&bpos,mesh->nsols,
1743 ncellSols,entities,type,size );
1744
1745 for (k=1; k<=mesh->nt; k++) {
1746 MMG5_pTria ptt = &mesh->tria[k];
1747 if ( !MG_EOK(ptt) ) continue;
1748
1749 for ( j=0; j<mesh->nsols; ++j ) {
1750 psl = *sol + j;
1751 if ( psl->entities==MMG5_Triangle ) {
1752 MMG2D_writeDoubleSol(psl,inm,bin,k,gmsh);
1753 }
1754 }
1755 fprintf(inm,"\n");
1756 }
1757
1758 MMG5_SAFE_FREE(type);
1759 MMG5_SAFE_FREE(size);
1760 MMG5_SAFE_FREE(entities);
1761
1762 /* End file */
1763 if ( !bin ) {
1764 fprintf(inm,"\n\nEnd\n");
1765 }
1766 else {
1767 binch = 54; //End
1768 fwrite(&binch,MMG5_SW,1,inm);
1769 }
1770 fclose(inm);
1771
1772 return 1;
1773}
1774
1775/* Custom version of Savemesh for debugging purpose */
1777 MMG5_pTria pt;
1778 MMG5_pEdge pa;
1779 MMG5_pPoint ppt,p0,p1,p2;
1780 MMG5_int k,np,nt,nc;
1781 FILE *out;
1782
1783 out = fopen(filename,"w");
1784
1785 np = nt = 0;
1786 /* Write Header */
1787 fprintf(out,"MeshVersionFormatted %d\n\nDimension %d\n\n",1,2);
1788
1789 /* Print vertices */
1790 for (k=1; k<=mesh->np; k++) {
1791 ppt = &mesh->point[k];
1792 if ( pack && MG_VOK(ppt) ) {
1793 np++;
1794 ppt->tmp = np;
1795 }
1796 else if ( !pack ) {
1797 np++;
1798 }
1799 }
1800
1801 fprintf(out,"Vertices\n %" MMG5_PRId "\n\n",np);
1802 for (k=1; k<=mesh->np; k++) {
1803 ppt = &mesh->point[k];
1804 if ( ( pack && MG_VOK(ppt) ) || !pack )
1805 fprintf(out,"%f %f %" MMG5_PRId "\n",ppt->c[0],ppt->c[1],ppt->ref);
1806 }
1807
1808 /* Print Triangles */
1809 for (k=1; k<=mesh->nt; k++) {
1810 pt = &mesh->tria[k];
1811 if ( MG_EOK(pt) ) nt++;
1812 }
1813
1814 fprintf(out,"Triangles\n %" MMG5_PRId "\n\n",nt);
1815 for (k=1; k<=mesh->nt; k++) {
1816 pt = &mesh->tria[k];
1817 if ( MG_EOK(pt) ) {
1818 p0 = &mesh->point[pt->v[0]];
1819 p1 = &mesh->point[pt->v[1]];
1820 p2 = &mesh->point[pt->v[2]];
1821 if ( pack ) {
1822 fprintf(out,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",p0->tmp,p1->tmp,p2->tmp,pt->ref);
1823 }
1824 else {
1825 fprintf(out,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",pt->v[0],pt->v[1],pt->v[2],pt->ref);
1826 }
1827 }
1828 }
1829
1830 /* Print Edges */
1831 if ( mesh->na ) {
1832 fprintf(out,"Edges\n %" MMG5_PRId "\n\n",mesh->na);
1833 for (k=1; k<=mesh->na; k++) {
1834 pa = &mesh->edge[k];
1835 p1 = &mesh->point[pa->a];
1836 p2 = &mesh->point[pa->b];
1837 if ( pack ) fprintf(out,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",p1->tmp,p2->tmp,pa->ref);
1838 else fprintf(out,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",pa->a,pa->b,pa->ref);
1839 }
1840 }
1841
1842 /* Print corners */
1843 nc = 0;
1844 for (k=1; k<=mesh->np; k++) {
1845 ppt = &mesh->point[k];
1846 if ( MG_VOK(ppt) && ppt->tag & MG_CRN ) nc++;
1847 }
1848
1849 if ( nc ) {
1850 fprintf(out,"Corners\n %" MMG5_PRId "\n\n",nc);
1851 for (k=1; k<=mesh->np; k++) {
1852 ppt = &mesh->point[k];
1853 if ( MG_VOK(ppt) && ppt->tag & MG_CRN ) {
1854 if ( pack ) fprintf(out,"%" MMG5_PRId "\n",ppt->tmp);
1855 else fprintf(out,"%" MMG5_PRId "\n",k);
1856 }
1857 }
1858 }
1859
1860 /* End keyword */
1861 fprintf(out,"End\n");
1862
1863 fclose(out);
1864
1865 return 1;
1866}
1867
1868/* Custom version of Savemet for debugging purpose */
1870 MMG5_pPoint ppt;
1871 MMG5_int k,np;
1872 char *ptr,typ=0,*data;
1873 FILE *out;
1874
1875 if ( met->size == 1 ) typ =1;
1876 else if ( met->size == 3 ) typ = 3;
1877
1878 MMG5_SAFE_CALLOC(data,strlen(filename)+6,char,return 0);
1879 strcpy(data,filename);
1880 ptr = strstr(data,".mesh");
1881 if ( ptr )
1882 *ptr = '\0';
1883
1884 strcat(data,".sol");
1885 out = fopen(data,"w");
1886
1887 MMG5_SAFE_FREE(data);
1888
1889 np = 0;
1890 for (k=1; k<=mesh->np; k++)
1891 mesh->point[k].tmp = 0;
1892
1893 /* Write Header */
1894 fprintf(out,"MeshVersionFormatted %d\n\nDimension %d\n\n",1,2);
1895
1896 /* Print vertices */
1897 for (k=1; k<=mesh->np; k++) {
1898 ppt = &mesh->point[k];
1899 if ( pack && MG_VOK(ppt) ) {
1900 np++;
1901 ppt->tmp = np;
1902 }
1903 else if ( !pack ) {
1904 np++;
1905 ppt->tmp = np;
1906 }
1907 }
1908
1909 fprintf(out,"SolAtVertices\n %" MMG5_PRId "\n%d %d\n\n",np,1,typ);
1910 for (k=1; k<=mesh->np; k++) {
1911 ppt = &mesh->point[k];
1912 if ( ( pack && MG_VOK(ppt) ) || !pack ) {
1913 if ( met->size == 1 )
1914 fprintf(out,"%f\n",met->m[k]);
1915 else if ( met->size == 3 )
1916 fprintf(out,"%f %f %f\n",met->m[3*k+0],met->m[3*k+1],met->m[3*k+2]);
1917 }
1918 }
1919
1920 /* End keyword */
1921 fprintf(out,"End\n");
1922
1923 fclose(out);
1924
1925 return 1;
1926}
1927
1928/* Save normal vector field for debugging purpose */
1930 MMG5_pPoint ppt;
1931 MMG5_int k,np;
1932 char *ptr,*data;
1933 FILE *out;
1934
1935 MMG5_SAFE_CALLOC(data,strlen(filename)+6,char,return 0);
1936 strcpy(data,filename);
1937 ptr = strstr(data,".mesh");
1938 if ( ptr )
1939 *ptr = '\0';
1940
1941 strcat(data,".nor.sol");
1942 out = fopen(data,"w");
1943
1944 MMG5_SAFE_FREE(data);
1945
1946 np = 0;
1947 for (k=1; k<=mesh->np; k++)
1948 mesh->point[k].tmp = 0;
1949
1950 /* Write Header */
1951 fprintf(out,"MeshVersionFormatted %d\n\nDimension %d\n\n",1,2);
1952
1953 /* Pack vertices or not for writing */
1954 for (k=1; k<=mesh->np; k++) {
1955 ppt = &mesh->point[k];
1956 if ( pack && MG_VOK(ppt) ) {
1957 np++;
1958 ppt->tmp = np;
1959 }
1960 else if ( !pack ) {
1961 np++;
1962 ppt->tmp = np;
1963 }
1964 }
1965
1966 fprintf(out,"SolAtVertices\n %" MMG5_PRId "\n%d %d\n\n",np,1,2);
1967 for (k=1; k<=mesh->np; k++) {
1968 ppt = &mesh->point[k];
1969 if ( ( pack && MG_VOK(ppt) ) || !pack ) {
1970 if ( MG_EDG(ppt->tag) && ! MG_SIN(ppt->tag) ) fprintf(out,"%f %f\n",ppt->n[0],ppt->n[1]);
1971 else fprintf(out,"%f %f\n",0.0,0.0);
1972 }
1973 }
1974
1975 /* End keyword */
1976 fprintf(out,"End\n");
1977
1978 fclose(out);
1979
1980 return 1;
1981}
1982
1983/* Save displacement field for debugging purpose */
1985 MMG5_pPoint ppt;
1986 MMG5_int k,np;
1987 char *ptr,*data;
1988 FILE *out;
1989
1990 MMG5_SAFE_CALLOC(data,strlen(filename)+6,char,return 0);
1991 strcpy(data,filename);
1992 ptr = strstr(data,".sol");
1993 if ( ptr )
1994 *ptr = '\0';
1995
1996 strcat(data,".disp.sol");
1997 out = fopen(data,"w");
1998 MMG5_SAFE_FREE(data);
1999
2000 np = 0;
2001 for (k=1; k<=mesh->np; k++)
2002 mesh->point[k].tmp = 0;
2003
2004 /* Write Header */
2005 fprintf(out,"MeshVersionFormatted %d\n\nDimension %d\n\n",1,2);
2006
2007 /* Pack vertices or not for writing */
2008 for (k=1; k<=mesh->np; k++) {
2009 ppt = &mesh->point[k];
2010 if ( pack && MG_VOK(ppt) ) {
2011 np++;
2012 ppt->tmp = np;
2013 }
2014 else if ( !pack ) {
2015 np++;
2016 ppt->tmp = np;
2017 }
2018 }
2019
2020 fprintf(out,"SolAtVertices\n %" MMG5_PRId "\n%d %d\n\n",np,1,2);
2021 for (k=1; k<=mesh->np; k++) {
2022 ppt = &mesh->point[k];
2023 if ( ( pack && MG_VOK(ppt) ) || !pack )
2024 fprintf(out,"%f %f\n",disp->m[2*(k-1)+1],disp->m[2*(k-1)+2]);
2025 }
2026
2027 /* End keyword */
2028 fprintf(out,"End\n");
2029
2030 fclose(out);
2031
2032 return 1;
2033}
2034
2035static inline
2037 FILE* inm;
2038 MMG5_pTria pt;
2039 MMG5_int k,i,ne;
2040 char *ptr,*data;
2041
2042 if ( !mesh->nt ) {
2043 return 1;
2044 }
2045
2046 if ( (!filename) || !(*filename) ) {
2048 }
2049 if ( (!filename) || !(*filename) ) {
2050 printf("\n ## Error: %s: unable to save a file without a valid filename\n.",
2051 __func__);
2052 return 0;
2053 }
2054
2055 /* Name of file */
2056 MMG5_SAFE_CALLOC(data,strlen(filename)+5,char,return 0);
2057 strcpy(data,filename);
2058 ptr = strstr(data,".node");
2059 if ( ptr ) {
2060 *ptr = '\0';
2061 }
2062
2063 /* Add .node ext */
2064 strcat(data,".ele");
2065 if( !(inm = fopen(data,"wb")) ) {
2066 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
2067 MMG5_SAFE_FREE(data);
2068 return 0;
2069 }
2070
2071 fprintf(stdout," %%%% %s OPENED\n",data);
2072 MMG5_SAFE_FREE(data);
2073
2074 ne = 0;
2075 for (k=1; k<=mesh->nt; k++) {
2076 pt = &mesh->tria[k];
2077 if ( !MG_EOK(pt) ) continue;
2078 ne++;
2079 }
2080
2081 /* Save elt number, node number per elt, 1 bdy marker */
2082 fprintf(inm, "%" MMG5_PRId " %d %d\n\n",ne,mesh->dim+1,1);
2083
2084 ne = 0;
2085 for ( k=1; k<=mesh->nt; ++k ) {
2086 pt = &mesh->tria[k];
2087 if ( MG_EOK(pt) ) {
2088 /* Save elt idx */
2089 fprintf(inm, "%" MMG5_PRId " ",++ne);
2090
2091 /* Save connectivity */
2092 for ( i=0; i<=mesh->dim; ++i ) {
2093 fprintf(inm, "%" MMG5_PRId " ",mesh->point[pt->v[i]].tmp);
2094 }
2095
2096 /* Save bdy marker */
2097 fprintf(inm, "%" MMG5_PRId "\n",pt->ref);
2098 }
2099 }
2100 fprintf(stdout," NUMBER OF ELEMENT %8" MMG5_PRId "\n",ne);
2101
2102 fclose(inm);
2103
2104 return 1;
2105}
2106
2107static inline
2109 FILE* inm;
2110 MMG5_pTria pt;
2111 MMG5_int k,i,ne,idx;
2112 char *ptr,*data;
2113
2114 if ( !mesh->nt ) {
2115 return 1;
2116 }
2117
2118 if ( (!filename) || !(*filename) ) {
2120 }
2121 if ( (!filename) || !(*filename) ) {
2122 printf("\n ## Error: %s: unable to save a file without a valid filename\n.",
2123 __func__);
2124 return 0;
2125 }
2126
2127 /* Name of file */
2128 MMG5_SAFE_CALLOC(data,strlen(filename)+7,char,return 0);
2129 strcpy(data,filename);
2130 ptr = strstr(data,".node");
2131 if ( ptr ) {
2132 *ptr = '\0';
2133 }
2134
2135 /* Add .node ext */
2136 strcat(data,".neigh");
2137 if( !(inm = fopen(data,"wb")) ) {
2138 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
2139 MMG5_SAFE_FREE(data);
2140 return 0;
2141 }
2142
2143 fprintf(stdout," %%%% %s OPENED\n",data);
2144 MMG5_SAFE_FREE(data);
2145
2146 if ( ! mesh->adja ) {
2147 if ( !MMG2D_hashTria(mesh) ) {
2148 printf("\n ## Error: %s: unable to compute triangle adjacencies\n.",__func__);
2149 return 0;
2150 }
2151 }
2152
2153 ne = 0;
2154 for (k=1; k<=mesh->nt; k++) {
2155 pt = &mesh->tria[k];
2156 if ( !MG_EOK(pt) ) continue;
2157 ne++;
2158 }
2159
2160 /* Save elt number, number of neighbors per elt */
2161 fprintf(inm, "%" MMG5_PRId " %d\n\n",ne,mesh->dim+1);
2162
2163 ne = 0;
2164 for ( k=1; k<=mesh->nt; ++k ) {
2165 pt = &mesh->tria[k];
2166 if ( MG_EOK(pt) ) {
2167 /* Save elt idx */
2168 fprintf(inm, "%" MMG5_PRId " ",++ne);
2169
2170 /* Save neighbors */
2171 for ( i=1; i<=mesh->dim+1; ++i ) {
2172 /* The triangle conventions is that no neighbors <=> -1 */
2173 idx = ( mesh->adja[3*(k-1)+i] > 0 ) ? mesh->adja[3*(k-1)+i]/3 : -1;
2174 fprintf(inm, "%" MMG5_PRId " ",idx);
2175 }
2176 /* Save bdy marker */
2177 fprintf(inm, "\n");
2178 }
2179 }
2180
2181 fclose(inm);
2182
2183 return 1;
2184}
2185
2186static inline
2188 MMG5_int nb_edges;
2189 int ier;
2190
2191 ier = MMG5_saveEdge(mesh,filename,".poly");
2192 if ( !ier ) {
2193 printf("\n ## Error: %s: unable to save boundary edges\n.",__func__);
2194 return 0;
2195 }
2196
2197 nb_edges = 0;
2198 ier = MMG2D_Get_numberOfNonBdyEdges( mesh, &nb_edges);
2199 if ( !ier ) {
2200 printf("\n ## Error: %s: unable to count and append internal edges\n.",__func__);
2201 return 0;
2202 }
2203
2204 ier = MMG5_saveEdge(mesh,filename,".edge");
2205 return ier;
2206}
2207
2208
2210
2211 if ( !MMG5_saveNode(mesh,filename) ) {
2212 return 0;
2213 }
2214
2215 if ( !MMG2D_saveEle(mesh,filename) ) {
2216 return 0;
2217 }
2218
2219 if ( !MMG2D_saveEdge(mesh,filename) ) {
2220 return 0;
2221 }
2222
2223 if ( !MMG2D_saveNeigh(mesh,filename) ) {
2224 return 0;
2225 }
2226
2227 return 1;
2228}
2229
2231 int ier=0;
2232 const char *filenameptr,*solnameptr;
2233 char *tmp,*soltmp;
2234
2235 if ( filename && strlen(filename) ) {
2236 filenameptr = filename;
2237 solnameptr = filename;
2238 }
2239 else if (mesh->namein && strlen(mesh->namein) ) {
2240 filenameptr = mesh->namein;
2241 if ( sol && strlen(sol->namein) ) {
2242 solnameptr = sol->namein;
2243 }
2244 else {
2245 solnameptr = mesh->namein;
2246 }
2247 }
2248 else {
2249 fprintf(stderr," ## Error: %s: please provide input file name"
2250 " (either in the mesh structure or as function argument).\n",
2251 __func__);
2252 return 0;
2253 }
2254
2255 MMG5_SAFE_MALLOC(tmp,strlen(filenameptr)+1,char,return 0);
2256 strcpy(tmp,filenameptr);
2257
2258 /* read mesh/sol files */
2259 char *ptr = MMG5_Get_filenameExt(tmp);
2260 int fmt = MMG5_Get_format(ptr,MMG5_FMT_MeditASCII);
2261
2262 int8_t savesolFile = 0;
2263
2264 switch ( fmt ) {
2265 case ( MMG5_FMT_GmshASCII ): case ( MMG5_FMT_GmshBinary ):
2267 break;
2268 case ( MMG5_FMT_VtkVtu ):
2270 break;
2271 case ( MMG5_FMT_VtkVtk ):
2273 break;
2274 case ( MMG5_FMT_Tetgen ):
2276 savesolFile = 1;
2277 break;
2278 default:
2280 savesolFile = 1;
2281 break;
2282 }
2283
2284 if ( ier && savesolFile ) {
2285 /* Medit or tetgen output: save the solution in a .sol file */
2286 if ( sol && sol->np ) {
2287 MMG5_SAFE_MALLOC(soltmp,strlen(solnameptr)+1,char,return 0);
2288 strcpy(soltmp,solnameptr);
2289
2290 if ( MMG2D_saveSol(mesh,sol,soltmp) == -1) {
2291 fprintf(stderr,"\n ## ERROR: WRONG DATA TYPE OR WRONG SOLUTION NUMBER.\n");
2292 ier = 0;
2293 }
2294 MMG5_SAFE_FREE(soltmp);
2295 }
2296 }
2297
2298 return ier;
2299}
const char * MMG5_Get_formatName(enum MMG5_Format fmt)
const char * MMG5_Get_entitiesName(enum MMG5_entities ent)
int MMG5_Get_format(char *ptr, int fmt)
char * MMG5_Get_filenameExt(char *filename)
int MMG2D_Set_solSize(MMG5_pMesh mesh, MMG5_pSol sol, int typEntity, MMG5_int np, int typSol)
int MMG2D_Set_inputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solin)
int ier
tmp[*strlen0]
MMG5_pMesh MMG5_pSol * sol
if(!ier) exit(EXIT_FAILURE)
MMG5_pMesh * mesh
MMG5_pMesh char * filename
int MMG2D_hashTria(MMG5_pMesh mesh)
Definition: hash_2d.c:35
int MMG5_chkMetricType(MMG5_pMesh mesh, int *type, int *entities, FILE *inm)
Definition: inout.c:2661
void MMG5_printSolStats(MMG5_pMesh mesh, MMG5_pSol *sol)
Definition: inout.c:2723
float MMG5_swapf(float sbin)
Definition: inout.c:80
int MMG5_saveSolAtTrianglesHeader(MMG5_pMesh mesh, FILE *inm, int ver, int bin, MMG5_int *bpos, int nsols, int nsolsAtTriangles, int *entities, int *type, int *size)
Definition: inout.c:2524
int MMG5_saveNode(MMG5_pMesh mesh, const char *filename)
Definition: inout.c:2741
int MMG5_saveEdge(MMG5_pMesh mesh, const char *filename, const char *ext)
Definition: inout.c:2816
int MMG5_saveSolHeader(MMG5_pMesh mesh, const char *filename, FILE **inm, int ver, int *bin, MMG5_int *bpos, MMG5_int np, int dim, int nsols, int *entities, int *type, int *size)
Definition: inout.c:2387
int MMG5_loadMshMesh_part1(MMG5_pMesh mesh, const char *filename, FILE **inm, long *posNodes, long *posElts, long **posNodeData, int *bin, int *iswp, MMG5_int *nelts, int *nsols)
Definition: inout.c:278
double MMG5_swapd(double sbin)
Definition: inout.c:97
int MMG5_loadMshMesh_part2(MMG5_pMesh mesh, MMG5_pSol *sol, FILE **inm, const long posNodes, const long posElts, const long *posNodeData, const int bin, const int iswp, const MMG5_int nelts, const int nsols)
Definition: inout.c:668
int MMG5_saveMshMesh(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename, int metricData)
Definition: inout.c:1587
void MMG5_printMetStats(MMG5_pMesh mesh, MMG5_pSol met)
Definition: inout.c:2705
int MMG5_loadSolHeader(const char *filename, int meshDim, FILE **inm, int *ver, int *bin, int *iswp, MMG5_int *np, int *dim, int *nsols, int **type, long *posnp, int imprim)
Definition: inout.c:2080
int MMG5_swapbin(int sbin)
Definition: inout.c:42
static int MMG2D_readDoubleSol(MMG5_pSol sol, FILE *inm, int bin, int iswp, MMG5_int pos)
Definition: inout_2d.c:873
int MMG2D_saveTetgenMesh(MMG5_pMesh mesh, const char *filename)
Definition: inout_2d.c:2209
int MMG2D_saveMesh(MMG5_pMesh mesh, const char *filename)
Definition: inout_2d.c:1096
int MMG2D_saveMshMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inout_2d.c:1535
int MMG2D_loadMshMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inout_2d.c:767
static int MMG2D_readFloatSol(MMG5_pSol sol, FILE *inm, int bin, int iswp, MMG5_int pos)
Definition: inout_2d.c:842
int MMG2D_saveMshMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inout_2d.c:1531
static int MMG2D_saveNeigh(MMG5_pMesh mesh, const char *filename)
Definition: inout_2d.c:2108
int MMG2D_savedisp_db(MMG5_pMesh mesh, MMG5_pSol disp, char *filename, int8_t pack)
Definition: inout_2d.c:1984
int MMG2D_loadGenericMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inout_2d.c:592
int MMG2D_saveGenericMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inout_2d.c:2230
int MMG2D_saveAllSols(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inout_2d.c:1673
int MMG2D_saveSol(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inout_2d.c:1612
int MMG2D_loadMesh(MMG5_pMesh mesh, const char *filename)
Definition: inout_2d.c:28
static int MMG2D_saveEle(MMG5_pMesh mesh, const char *filename)
Definition: inout_2d.c:2036
int MMG2D_savemesh_db(MMG5_pMesh mesh, char *filename, int8_t pack)
Definition: inout_2d.c:1776
int MMG2D_loadSol(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inout_2d.c:900
int MMG2D_loadMshMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inout_2d.c:701
int MMG2D_savemet_db(MMG5_pMesh mesh, MMG5_pSol met, char *filename, int8_t pack)
Definition: inout_2d.c:1869
int MMG2D_2dMeshCheck(MMG5_pMesh mesh)
Definition: inout_2d.c:674
static int MMG2D_saveEdge(MMG5_pMesh mesh, const char *filename)
Definition: inout_2d.c:2187
int MMG2D_savenor_db(MMG5_pMesh mesh, char *filename, int8_t pack)
Definition: inout_2d.c:1929
int MMG2D_loadAllSols(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inout_2d.c:991
static void MMG2D_writeDoubleSol(MMG5_pSol sol, FILE *inm, int bin, MMG5_int pos, int gmsh)
Definition: inout_2d.c:1550
int MMG2D_loadVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
int MMG2D_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
int MMG2D_saveVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
int MMG2D_loadVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
API headers for the mmg2d library.
LIBMMG2D_EXPORT int MMG2D_Get_numberOfNonBdyEdges(MMG5_pMesh mesh, MMG5_int *nb_edges)
int MMG2D_zaldy(MMG5_pMesh mesh)
Definition: zaldy_2d.c:303
@ MMG5_FMT_MeditBinary
Definition: libmmgtypes.h:236
@ MMG5_FMT_Tetgen
Definition: libmmgtypes.h:244
@ MMG5_FMT_VtkVtk
Definition: libmmgtypes.h:243
@ MMG5_FMT_GmshBinary
Definition: libmmgtypes.h:238
@ MMG5_FMT_GmshASCII
Definition: libmmgtypes.h:237
@ MMG5_FMT_MeditASCII
Definition: libmmgtypes.h:235
@ MMG5_FMT_VtkVtu
Definition: libmmgtypes.h:241
#define MMG5_NSOLS_MAX
Definition: libmmgtypes.h:181
@ MMG5_Noentity
Definition: libmmgtypes.h:223
@ MMG5_Vertex
Definition: libmmgtypes.h:224
@ MMG5_Triangle
Definition: libmmgtypes.h:226
#define MG_REQ
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
#define MMG5_EPSOK
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_SD
#define MMG5_ADD_MEM(mesh, size, message, law)
#define MG_EDG(tag)
#define MG_SIN(tag)
#define MMG_FREAD(ptr, size, count, stream)
#define MMG5_SW
#define MG_VOK(ppt)
#define MG_CRN
#define MG_BDY
#define MMG5_SAFE_MALLOC(ptr, size, type, law)
#define MG_NOSURF
#define MMG5_FILESTR_LGTH
#define MG_REF
#define MMG5_SAFE_FREE(ptr)
#define MMG_FSCANF(stream, format,...)
#define MMG5_DEL_MEM(mesh, ptr)
Structure to store edges of a MMG mesh.
Definition: libmmgtypes.h:305
MMG5_int b
Definition: libmmgtypes.h:306
MMG5_int ref
Definition: libmmgtypes.h:307
int16_t tag
Definition: libmmgtypes.h:310
MMG5_int a
Definition: libmmgtypes.h:306
uint8_t nosurf
Definition: libmmgtypes.h:546
MMG mesh structure.
Definition: libmmgtypes.h:605
MMG5_int ntmax
Definition: libmmgtypes.h:612
MMG5_pQuad quadra
Definition: libmmgtypes.h:648
MMG5_Info info
Definition: libmmgtypes.h:651
MMG5_int ne
Definition: libmmgtypes.h:612
MMG5_pPoint point
Definition: libmmgtypes.h:641
char * nameout
Definition: libmmgtypes.h:653
MMG5_int * adja
Definition: libmmgtypes.h:624
MMG5_int nquad
Definition: libmmgtypes.h:613
MMG5_int npmax
Definition: libmmgtypes.h:612
MMG5_int xp
Definition: libmmgtypes.h:620
MMG5_int nt
Definition: libmmgtypes.h:612
MMG5_pTria tria
Definition: libmmgtypes.h:647
MMG5_int np
Definition: libmmgtypes.h:612
MMG5_pEdge edge
Definition: libmmgtypes.h:649
MMG5_int nti
Definition: libmmgtypes.h:612
MMG5_int npi
Definition: libmmgtypes.h:612
MMG5_int nai
Definition: libmmgtypes.h:612
MMG5_int nprism
Definition: libmmgtypes.h:613
MMG5_int na
Definition: libmmgtypes.h:612
char * namein
Definition: libmmgtypes.h:652
Structure to store points of a MMG mesh.
Definition: libmmgtypes.h:270
double n[3]
Definition: libmmgtypes.h:272
int16_t tag
Definition: libmmgtypes.h:284
MMG5_int tmp
Definition: libmmgtypes.h:280
double c[3]
Definition: libmmgtypes.h:271
MMG5_int ref
Definition: libmmgtypes.h:278
MMG5_int ref
Definition: libmmgtypes.h:368
int16_t tag[4]
Definition: libmmgtypes.h:372
MMG5_int v[4]
Definition: libmmgtypes.h:367
int entities
Definition: libmmgtypes.h:670
double * m
Definition: libmmgtypes.h:671
MMG5_int edg[3]
Definition: libmmgtypes.h:339
int16_t tag[3]
Definition: libmmgtypes.h:342
MMG5_int ref
Definition: libmmgtypes.h:335
MMG5_int v[3]
Definition: libmmgtypes.h:334