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