Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
inout_3d.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
36#include "libmmg3d.h"
37#include "libmmg3d_private.h"
38
53int MMG3D_openMesh(int imprim,const char *filename,FILE **inm,int *bin,
54 char *modeASCII, char* modeBIN) {
55 char *ptr,*data;
56 int out,ier;
57
58 out = (strchr(modeASCII,'w') != NULL) ? 1 : 0;
59 ier = out ? 0 : -1;
60
61 *bin = 0;
62 MMG5_SAFE_CALLOC(data,strlen(filename)+7,char,return ier);
63
64 strcpy(data,filename);
65 ptr = strstr(data,".mesh");
66
67 if ( !ptr ) {
68 /* data contains the filename without extension */
69 strcat(data,".meshb");
70 if( !(*inm = fopen(data,modeBIN)) ) {
71 /* our file is not a .meshb file, try with .mesh ext */
72 ptr = strstr(data,".mesh");
73 *ptr = '\0';
74 strcat(data,".mesh");
75 if( !(*inm = fopen(data,modeASCII)) ) {
76 MMG5_SAFE_FREE(data);
77 return 0;
78 }
79 }
80 else *bin = 1;
81 }
82 else {
83 ptr = strstr(data,".meshb");
84 if ( ptr ) {
85 *bin = 1;
86 if( !(*inm = fopen(data,modeBIN)) ) {
87 if ( out ) {
88 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
89 }
90 MMG5_SAFE_FREE(data);
91 return 0;
92 }
93 }
94 else {
95 if( !(*inm = fopen(data,modeASCII)) ) {
96 if ( out ) {
97 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
98 }
99 MMG5_SAFE_FREE(data);
100 return 0;
101 }
102 }
103 }
104
105 if ( imprim >= 0 ) {
106 fprintf(stdout," %%%% %s OPENED\n",data);
107 }
108 MMG5_SAFE_FREE(data);
109
110 return 1;
111}
112
113int MMG3D_loadMesh_opened(MMG5_pMesh mesh,FILE *inm,int bin) {
114 MMG5_pTetra pt;
115 MMG5_pPrism pp;
116 MMG5_pTria pt1;
117 MMG5_pQuad pq1;
118 MMG5_pEdge pa;
119 MMG5_pPoint ppt;
120 double *norm,*n,dd;
121 float fc;
122 long posnp,posnt,posne,posned,posncor,posnpreq,posntreq,posnereq,posnedreq;
123 long posnr,posnprism,posnormal,posnc1,posnq,posnqreq;
124 long posnppar,posntpar,posnepar,posnedpar,posnqpar;
125 MMG5_int npreq,ntreq,nereq,nedreq,nqreq,ncor,ned,ng;
126 MMG5_int nppar,nedpar,ntpar,nqpar,nepar;
127 MMG5_int k,ip,idn;
128 int i,bdim,binch,iswp,bpos;
129 MMG5_int na,nr,ia,aux,nref,ref;
130 char chaine[MMG5_FILESTR_LGTH],strskip[MMG5_FILESTR_LGTH];
131
132 posnp = posnt = posne = posncor = 0;
133 posnpreq = posntreq = posnereq = posnqreq = posned = posnedreq = posnr = 0;
134 posnprism = posnormal= posnc1 = posnq = 0;
135 posnppar = posnedpar = posntpar = posnqpar = posnepar = 0;
136 ncor = ned = npreq = ntreq = nqreq = nereq = nedreq = nr = ng = 0;
137 nppar = nedpar = ntpar = nqpar = nepar = 0;
138 iswp = 0;
139 bpos = ia = idn = ip = 0;
140 mesh->np = mesh->nt = mesh->ne = 0;
141 nref = 0;
142
143
144 if (!bin) {
145 strcpy(chaine,"D");
146 while(fscanf(inm,"%127s",&chaine[0])!=EOF && strncmp(chaine,"End",strlen("End")) ) {
147 if ( chaine[0] == '#' ) {
148 while(1){ // skip until end of line or file
149 char *s = fgets(strskip,MMG5_FILESTR_LGTH,inm);
150 if(!s) break; // nothing could be read
151 if(s[strlen(s)-1]=='\n') break; // end of line
152 }
153 continue;
154 }
155
156 if(!strncmp(chaine,"MeshVersionFormatted",strlen("MeshVersionFormatted"))) {
157 MMG_FSCANF(inm,"%d",&mesh->ver);
158 continue;
159 } else if(!strncmp(chaine,"Dimension",strlen("Dimension"))) {
160 MMG_FSCANF(inm,"%d",&mesh->dim);
161 if(mesh->dim!=3) {
162 fprintf(stderr,"BAD DIMENSION : %d\n",mesh->dim);
163 return -1;
164 }
165 continue;
166 } else if(!strncmp(chaine,"Vertices",strlen("Vertices"))) {
167 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->npi);
168 posnp = ftell(inm);
169 continue;
170 } else if(!strncmp(chaine,"RequiredVertices",strlen("RequiredVertices"))) {
171 MMG_FSCANF(inm,"%" MMG5_PRId "",&npreq);
172 posnpreq = ftell(inm);
173 continue;
174 } else if(!strncmp(chaine,"ParallelVertices",strlen("ParallelVertices"))) {
175 MMG_FSCANF(inm,"%" MMG5_PRId "",&nppar);
176 posnppar = ftell(inm);
177 continue;
178 } else if(!strncmp(chaine,"Triangles",strlen("Triangles"))) {
179 if ( !strncmp(chaine,"TrianglesP",strlen("TrianglesP")) ) continue;
180 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->nti);
181 posnt = ftell(inm);
182 continue;
183 } else if(!strncmp(chaine,"RequiredTriangles",strlen("RequiredTriangles"))) {
184 MMG_FSCANF(inm,"%" MMG5_PRId "",&ntreq);
185 posntreq = ftell(inm);
186 continue;
187 } else if(!strncmp(chaine,"ParallelTriangles",strlen("ParallelTriangles"))) {
188 MMG_FSCANF(inm,"%" MMG5_PRId "",&ntpar);
189 posntpar = ftell(inm);
190 continue;
191 }
192 else if(!strncmp(chaine,"Quadrilaterals",strlen("Quadrilaterals"))) {
193 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->nquad);
194 posnq = ftell(inm);
195 continue;
196 } else if(!strncmp(chaine,"RequiredQuadrilaterals",strlen("RequiredQuadrilaterals"))) {
197 MMG_FSCANF(inm,"%" MMG5_PRId "",&nqreq);
198 posnqreq = ftell(inm);
199 continue;
200 } else if(!strncmp(chaine,"ParallelQuadrilaterals",strlen("ParallelQuadrilaterals"))) {
201 MMG_FSCANF(inm,"%" MMG5_PRId "",&nqpar);
202 posnqpar = ftell(inm);
203 continue;
204 } else if(!strncmp(chaine,"Tetrahedra",strlen("Tetrahedra"))) {
205 if ( !strncmp(chaine,"TetrahedraP",strlen("TetrahedraP")) ) continue;
206 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->nei);
207 posne = ftell(inm);
208 continue;
209 } else if((!strncmp(chaine,"Prisms",strlen("Prisms")))||
210 (!strncmp(chaine,"Pentahedra",strlen("Pentahedra")))) {
211 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->nprism);
212 posnprism = ftell(inm);
213 continue;
214 } else if(!strncmp(chaine,"RequiredTetrahedra",strlen("RequiredTetrahedra"))) {
215 MMG_FSCANF(inm,"%" MMG5_PRId "",&nereq);
216 posnereq = ftell(inm);
217 continue;
218 } else if(!strncmp(chaine,"ParallelTetrahedra",strlen("ParallelTetrahedra"))) {
219 MMG_FSCANF(inm,"%" MMG5_PRId "",&nepar);
220 posnepar = ftell(inm);
221 continue;
222 } else if(!strncmp(chaine,"Corners",strlen("Corners"))) {
223 MMG_FSCANF(inm,"%" MMG5_PRId "",&ncor);
224 posncor = ftell(inm);
225 continue;
226 } else if(!strncmp(chaine,"Edges",strlen("Edges"))) {
227 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->nai);
228 posned = ftell(inm);
229 continue;
230 } else if(!strncmp(chaine,"RequiredEdges",strlen("RequiredEdges"))) {
231 MMG_FSCANF(inm,"%" MMG5_PRId "",&nedreq);
232 posnedreq = ftell(inm);
233 continue;
234 } else if(!strncmp(chaine,"ParallelEdges",strlen("ParallelEdges"))) {
235 MMG_FSCANF(inm,"%" MMG5_PRId "",&nedpar);
236 posnedpar = ftell(inm);
237 continue;
238 } else if(!strncmp(chaine,"Ridges",strlen("Ridges"))) {
239 MMG_FSCANF(inm,"%" MMG5_PRId "",&nr);
240 posnr = ftell(inm);
241 continue;
242 } else if(!ng && !strncmp(chaine,"Normals",strlen("Normals"))) {
243 MMG_FSCANF(inm,"%" MMG5_PRId "",&ng);
244 posnormal = ftell(inm);
245 continue;
246 } else if(!strncmp(chaine,"NormalAtVertices",strlen("NormalAtVertices"))) {
247 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->nc1);
248 posnc1 = ftell(inm);
249 continue;
250 }
251 }
252 } else { //binary file
253 bdim = 0;
254 MMG_FREAD(&mesh->ver,MMG5_SW,1,inm);
255 iswp=0;
256 if(mesh->ver==16777216)
257 iswp=1;
258 else if(mesh->ver!=1) {
259 fprintf(stderr,"BAD FILE ENCODING\n");
260 }
261 MMG_FREAD(&mesh->ver,MMG5_SW,1,inm);
262 if(iswp) mesh->ver = MMG5_swapbin(mesh->ver);
263 while(fread(&binch,MMG5_SW,1,inm)!=0 && binch!=54 ) {
264 if(iswp) binch=MMG5_swapbin(binch);
265 if(binch==54) break;
266 if(!bdim && binch==3) { //Dimension
267 MMG_FREAD(&bdim,MMG5_SW,1,inm); //NulPos=>20
268 if(iswp) bdim=MMG5_swapbin(bdim);
269 MMG_FREAD(&bdim,MMG5_SW,1,inm);
270 if(iswp) bdim=MMG5_swapbin(bdim);
271 mesh->dim = bdim;
272 if(bdim!=3) {
273 fprintf(stderr,"BAD MESH DIMENSION : %d\n",mesh->dim);
274 fprintf(stderr," Exit program.\n");
275 return -1;
276 }
277 continue;
278 } else if(!mesh->npi && binch==4) { //Vertices
279 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
280 if(iswp) bpos=MMG5_swapbin(bpos);
281 MMG_FREAD(&mesh->npi,MMG5_SW,1,inm);
282 if(iswp) mesh->npi=MMG5_swapbin(mesh->npi);
283 posnp = ftell(inm);
284 rewind(inm);
285 fseek(inm,bpos,SEEK_SET);
286 continue;
287 } else if(binch==15) { //RequiredVertices
288 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
289 if(iswp) bpos=MMG5_swapbin(bpos);
290 MMG_FREAD(&npreq,MMG5_SW,1,inm);
291 if(iswp) npreq=MMG5_swapbin(npreq);
292 posnpreq = ftell(inm);
293 rewind(inm);
294 fseek(inm,bpos,SEEK_SET);
295 continue;
296 } else if(!mesh->nti && binch==6) {//Triangles
297 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
298 if(iswp) bpos=MMG5_swapbin(bpos);
299 MMG_FREAD(&mesh->nti,MMG5_SW,1,inm);
300 if(iswp) mesh->nti=MMG5_swapbin(mesh->nti);
301 posnt = ftell(inm);
302 rewind(inm);
303 fseek(inm,bpos,SEEK_SET);
304 continue;
305 } else if(binch==17) { //RequiredTriangles
306 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
307 if(iswp) bpos=MMG5_swapbin(bpos);
308 MMG_FREAD(&ntreq,MMG5_SW,1,inm);
309 if(iswp) ntreq=MMG5_swapbin(ntreq);
310 posntreq = ftell(inm);
311 rewind(inm);
312 fseek(inm,bpos,SEEK_SET);
313 continue;
314 }
315 else if(!mesh->nquad && binch==7) {//Quadrilaterals
316 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
317 if(iswp) bpos=MMG5_swapbin(bpos);
318 MMG_FREAD(&mesh->nquad,MMG5_SW,1,inm);
319 if(iswp) mesh->nquad=MMG5_swapbin(mesh->nquad);
320 posnq = ftell(inm);
321 rewind(inm);
322 fseek(inm,bpos,SEEK_SET);
323 continue;
324 } else if(binch==18) { //RequiredQuadrilaterals
325 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
326 if(iswp) bpos=MMG5_swapbin(bpos);
327 MMG_FREAD(&nqreq,MMG5_SW,1,inm);
328 if(iswp) nqreq=MMG5_swapbin(nqreq);
329 posnqreq = ftell(inm);
330 rewind(inm);
331 fseek(inm,bpos,SEEK_SET);
332 continue;
333 } else if(!mesh->nei && binch==8) {//Tetra
334 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
335 if(iswp) bpos=MMG5_swapbin(bpos);
336 MMG_FREAD(&mesh->nei,MMG5_SW,1,inm);
337 if(iswp) mesh->nei=MMG5_swapbin(mesh->nei);
338 posne = ftell(inm);
339 rewind(inm);
340 fseek(inm,bpos,SEEK_SET);
341 continue;
342 } else if(!mesh->nprism && binch==9) {//Prism
343 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
344 if(iswp) bpos=MMG5_swapbin(bpos);
345 MMG_FREAD(&mesh->nprism,MMG5_SW,1,inm);
346 if(iswp) mesh->nprism=MMG5_swapbin(mesh->nprism);
347 posnprism = ftell(inm);
348 rewind(inm);
349 fseek(inm,bpos,SEEK_SET);
350 continue;
351 } else if(binch==12) { //RequiredTetra
352 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
353 if(iswp) bpos=MMG5_swapbin(bpos);
354 MMG_FREAD(&nereq,MMG5_SW,1,inm);
355 if(iswp) nereq=MMG5_swapbin(nereq);
356 posnereq = ftell(inm);
357 rewind(inm);
358 fseek(inm,bpos,SEEK_SET);
359 continue;
360 } else if(!ncor && binch==13) { //Corners
361 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
362 if(iswp) bpos=MMG5_swapbin(bpos);
363 MMG_FREAD(&ncor,MMG5_SW,1,inm);
364 if(iswp) ncor=MMG5_swapbin(ncor);
365 posncor = ftell(inm);
366 rewind(inm);
367 fseek(inm,bpos,SEEK_SET);
368 continue;
369 } else if(!mesh->nai && binch==5) { //Edges
370 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
371 if(iswp) bpos=MMG5_swapbin(bpos);
372 MMG_FREAD(&mesh->nai,MMG5_SW,1,inm);
373 if(iswp) mesh->nai=MMG5_swapbin(mesh->nai);
374 posned = ftell(inm);
375 rewind(inm);
376 fseek(inm,bpos,SEEK_SET);
377 continue;
378 } else if(binch==16) { //RequiredEdges
379 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
380 if(iswp) bpos=MMG5_swapbin(bpos);
381 MMG_FREAD(&nedreq,MMG5_SW,1,inm);
382 if(iswp) nedreq=MMG5_swapbin(nedreq);
383 posnedreq = ftell(inm);
384 rewind(inm);
385 fseek(inm,bpos,SEEK_SET);
386 continue;
387 } else if(binch==14) { //Ridges
388 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
389 if(iswp) bpos=MMG5_swapbin(bpos);
390 MMG_FREAD(&nr,MMG5_SW,1,inm);
391 if(iswp) nr=MMG5_swapbin(nr);
392 posnr = ftell(inm);
393 rewind(inm);
394 fseek(inm,bpos,SEEK_SET);
395 continue;
396 } else if(!ng && binch==60) { //Normals
397 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
398 if(iswp) bpos=MMG5_swapbin(bpos);
399 MMG_FREAD(&ng,MMG5_SW,1,inm);
400 if(iswp) ng=MMG5_swapbin(ng);
401 posnormal = ftell(inm);
402 rewind(inm);
403 fseek(inm,bpos,SEEK_SET);
404 continue;
405 } else if(binch==20) { //NormalAtVertices
406 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
407 if(iswp) bpos=MMG5_swapbin(bpos);
408 MMG_FREAD(&mesh->nc1,MMG5_SW,1,inm);
409 if(iswp) mesh->nc1=MMG5_swapbin(mesh->nc1);
410 posnc1 = ftell(inm);
411 rewind(inm);
412 fseek(inm,bpos,SEEK_SET);
413 continue;
414 } else {
415 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
416 if(iswp) bpos=MMG5_swapbin(bpos);
417
418 rewind(inm);
419 fseek(inm,bpos,SEEK_SET);
420 }
421 }
422 }
423
424 if ( !mesh->npi || !mesh->nei ) {
425 fprintf(stderr," ** MISSING DATA.\n");
426 fprintf(stderr," Check that your mesh contains points and tetrahedra.\n");
427 fprintf(stderr," Exit program.\n");
428 return -1;
429 }
430 /* memory allocation */
431 mesh->np = mesh->npi;
432 mesh->nt = mesh->nti;
433 mesh->ne = mesh->nei;
434 mesh->na = mesh->nai;
435 if ( !MMG3D_zaldy(mesh) ) return 0;
436 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt || mesh->nemax < mesh->ne) {
437 return -1;
438 }
439
440 rewind(inm);
441 fseek(inm,posnp,SEEK_SET);
442 for (k=1; k<=mesh->np; k++) {
443 ppt = &mesh->point[k];
444 if (mesh->ver < 2) { /*float*/
445 if (!bin) {
446 for (i=0 ; i<3 ; i++) {
447 MMG_FSCANF(inm,"%f",&fc);
448 ppt->c[i] = (double) fc;
449 }
450 MMG_FSCANF(inm,"%" MMG5_PRId "",&ppt->ref);
451 } else {
452 for (i=0 ; i<3 ; i++) {
453 MMG_FREAD(&fc,MMG5_SW,1,inm);
454 if(iswp) fc=MMG5_swapf(fc);
455 ppt->c[i] = (double) fc;
456 }
457 MMG_FREAD(&ppt->ref,MMG5_SW,1,inm);
458 if(iswp) ppt->ref=MMG5_swapbin(ppt->ref);
459 }
460 } else {
461 if (!bin) {
462 MMG_FSCANF(inm,"%lf %lf %lf %" MMG5_PRId "",&ppt->c[0],&ppt->c[1],&ppt->c[2],&ppt->ref);
463 }
464 else {
465 for (i=0 ; i<3 ; i++) {
466 MMG_FREAD(&ppt->c[i],MMG5_SD,1,inm);
467 if(iswp) ppt->c[i]=MMG5_swapd(ppt->c[i]);
468 }
469 MMG_FREAD(&ppt->ref,MMG5_SW,1,inm);
470 if(iswp) ppt->ref=MMG5_swapbin(ppt->ref);
471 }
472 }
473
474 if ( ppt->ref < 0 ) {
475 ppt->ref = -ppt->ref;
476 ++nref;
477 }
478
479 ppt->tag = MG_NUL;
480 ppt->tmp = 0;
481 }
482 /* get required vertices */
483 if(npreq) {
484 rewind(inm);
485 fseek(inm,posnpreq,SEEK_SET);
486 for (k=1; k<=npreq; k++) {
487 if(!bin) {
488 MMG_FSCANF(inm,"%d",&i);
489 }
490 else {
491 MMG_FREAD(&i,MMG5_SW,1,inm);
492 if(iswp) i=MMG5_swapbin(i);
493 }
494 if(i>mesh->np) {
495 fprintf(stderr,"\n ## Warning: %s: required Vertices number %8d"
496 " ignored.\n",__func__,i);
497 } else {
498 ppt = &mesh->point[i];
499 ppt->tag |= MG_REQ;
500 ppt->tag &= ~MG_NUL;
501 }
502 }
503 }
504 /* get parallel vertices */
505 if(nppar) {
506 rewind(inm);
507 fseek(inm,posnppar,SEEK_SET);
508 for (k=1; k<=nppar; k++) {
509 if(!bin) {
510 MMG_FSCANF(inm,"%d",&i);
511 }
512 else {
513 MMG_FREAD(&i,MMG5_SW,1,inm);
514 if(iswp) i=MMG5_swapbin(i);
515 }
516 if(i>mesh->np) {
517 fprintf(stderr,"\n ## Warning: %s: parallel Vertices number %8d"
518 " ignored.\n",__func__,i);
519 } else {
520 ppt = &mesh->point[i];
521 ppt->tag |= MG_PARBDY;
522 }
523 }
524 }
525
526
527 /* get corners */
528 if(ncor) {
529 rewind(inm);
530 fseek(inm,posncor,SEEK_SET);
531 for (k=1; k<=ncor; k++) {
532 if(!bin) {
533 MMG_FSCANF(inm,"%d",&i);
534 }
535 else {
536 MMG_FREAD(&i,MMG5_SW,1,inm);
537 if(iswp) i=MMG5_swapbin(i);
538 }
539 if(i>mesh->np) {
540 fprintf(stderr,"\n ## Warning: %s: corner number %8d ignored.\n",
541 __func__,i);
542 } else {
543 ppt = &mesh->point[i];
544 ppt->tag |= MG_CRN;
545 }
546 }
547 }
548
549 /* read mesh triangles */
550 if ( mesh->nt ) {
551 rewind(inm);
552 fseek(inm,posnt,SEEK_SET);
553 /* Skip triangles with mesh->info.isoref refs */
554 for (k=1; k<=mesh->nt; k++) {
555 pt1 = &mesh->tria[k];
556 if (!bin) {
557 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&pt1->v[0],&pt1->v[1],&pt1->v[2],&pt1->ref);
558 }
559 else {
560 for (i=0 ; i<3 ; i++) {
561 MMG_FREAD(&pt1->v[i],MMG5_SW,1,inm);
562 if(iswp) pt1->v[i]=MMG5_swapbin(pt1->v[i]);
563 }
564 MMG_FREAD(&pt1->ref,MMG5_SW,1,inm);
565 if(iswp) pt1->ref=MMG5_swapbin(pt1->ref);
566 }
567 if ( pt1->ref < 0 ) {
568 pt1->ref = -pt1->ref;
569 ++nref;
570 }
571 }
572 /* get required triangles */
573 if(ntreq) {
574 rewind(inm);
575 fseek(inm,posntreq,SEEK_SET);
576 for (k=1; k<=ntreq; k++) {
577 if(!bin) {
578 MMG_FSCANF(inm,"%d",&i);
579 }
580 else {
581 MMG_FREAD(&i,MMG5_SW,1,inm);
582 if(iswp) i=MMG5_swapbin(i);
583 }
584 if ( i>mesh->nt ) {
585 fprintf(stderr,"\n ## Warning: %s: required triangle number %8d"
586 " ignored.\n",__func__,i);
587 } else {
588 pt1 = &mesh->tria[i];
589 pt1->tag[0] |= MG_REQ;
590 pt1->tag[1] |= MG_REQ;
591 pt1->tag[2] |= MG_REQ;
592 }
593 }
594 }
595 /* get parallel triangles */
596 if(ntpar) {
597 rewind(inm);
598 fseek(inm,posntpar,SEEK_SET);
599 for (k=1; k<=ntpar; k++) {
600 if(!bin) {
601 MMG_FSCANF(inm,"%d",&i);
602 }
603 else {
604 MMG_FREAD(&i,MMG5_SW,1,inm);
605 if(iswp) i=MMG5_swapbin(i);
606 }
607 if ( i>mesh->nt ) {
608 fprintf(stderr,"\n ## Warning: %s: parallel triangle number %8d"
609 " ignored.\n",__func__,i);
610 } else {
611 pt1 = &mesh->tria[i];
612 pt1->tag[0] |= MG_PARBDY;
613 pt1->tag[1] |= MG_PARBDY;
614 pt1->tag[2] |= MG_PARBDY;
615 }
616 }
617 }
618 } //end if mesh->nt
619
620 /* read mesh quadrilaterals */
621 if ( mesh->nquad ) {
622 rewind(inm);
623 fseek(inm,posnq,SEEK_SET);
624
625 for (k=1; k<=mesh->nquad; k++) {
626 pq1 = &mesh->quadra[k];
627 if (!bin) {
628 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&pq1->v[0],&pq1->v[1],&pq1->v[2],
629 &pq1->v[3],&pq1->ref);
630 }
631 else {
632 for (i=0 ; i<4 ; i++) {
633 MMG_FREAD(&pq1->v[i],MMG5_SW,1,inm);
634 if(iswp) pq1->v[i]=MMG5_swapbin(pq1->v[i]);
635 }
636 MMG_FREAD(&pq1->ref,MMG5_SW,1,inm);
637 if(iswp) pq1->ref=MMG5_swapbin(pq1->ref);
638 }
639 if ( pq1->ref < 0 ) {
640 pq1->ref = -pq1->ref;
641 ++nref;
642 }
643 }
644
645 /* get required quadrilaterals */
646 if(nqreq) {
647 rewind(inm);
648 fseek(inm,posnqreq,SEEK_SET);
649 for (k=1; k<=nqreq; k++) {
650 if(!bin) {
651 MMG_FSCANF(inm,"%d",&i);
652 }
653 else {
654 MMG_FREAD(&i,MMG5_SW,1,inm);
655 if(iswp) i=MMG5_swapbin(i);
656 }
657 if ( i>mesh->nquad ) {
658 fprintf(stderr,"\n ## Warning: %s: required quadrilaterals number"
659 " %8d ignored.\n",__func__,i);
660 } else {
661 pq1 = &mesh->quadra[i];
662 pq1->tag[0] |= MG_REQ;
663 pq1->tag[1] |= MG_REQ;
664 pq1->tag[2] |= MG_REQ;
665 pq1->tag[3] |= MG_REQ;
666 }
667 }
668 }
669
670 /* get parallel quadrilaterals */
671 if(nqpar) {
672 rewind(inm);
673 fseek(inm,posnqpar,SEEK_SET);
674 for (k=1; k<=nqpar; k++) {
675 if(!bin) {
676 MMG_FSCANF(inm,"%d",&i);
677 }
678 else {
679 MMG_FREAD(&i,MMG5_SW,1,inm);
680 if(iswp) i=MMG5_swapbin(i);
681 }
682 if ( i>mesh->nquad ) {
683 fprintf(stderr,"\n ## Warning: %s: parallel quadrilaterals number"
684 " %8d ignored.\n",__func__,i);
685 } else {
686 pq1 = &mesh->quadra[i];
687 pq1->tag[0] |= MG_PARBDY;
688 pq1->tag[1] |= MG_PARBDY;
689 pq1->tag[2] |= MG_PARBDY;
690 pq1->tag[3] |= MG_PARBDY;
691 }
692 }
693 }
694 } //end if mesh->nquad
695
696 /* read mesh edges */
697 if ( mesh->na ) {
698 na = mesh->na;
699
700 rewind(inm);
701 fseek(inm,posned,SEEK_SET);
702
703 for (k=1; k<=na; k++) {
704 pa = &mesh->edge[k];
705 if (!bin) {
706 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&pa->a,&pa->b,&pa->ref);
707 }
708 else {
709 MMG_FREAD(&pa->a,MMG5_SW,1,inm);
710 if(iswp) pa->a=MMG5_swapbin(pa->a);
711 MMG_FREAD(&pa->b,MMG5_SW,1,inm);
712 if(iswp) pa->b=MMG5_swapbin(pa->b);
713 MMG_FREAD(&pa->ref,MMG5_SW,1,inm);
714 if(iswp) pa->ref=MMG5_swapbin(pa->ref);
715 }
716 pa->tag |= MG_REF;
717 if ( pa->ref < 0 ) {
718 pa->ref = -pa->ref;
719 ++nref;
720 }
721 }
722
723 /* get ridges */
724 if ( nr ) {
725 rewind(inm);
726 fseek(inm,posnr,SEEK_SET);
727 for (k=1; k<=nr; k++) {
728 if(!bin) {
729 MMG_FSCANF(inm,"%" MMG5_PRId "",&ia);
730 }
731 else {
732 MMG_FREAD(&ia,MMG5_SW,1,inm);
733 if(iswp) ia=MMG5_swapbin(ia);
734 }
735 if(ia>na) {
736 fprintf(stderr,"\n ## Warning: %s: ridge number %8" MMG5_PRId " ignored.\n",
737 __func__,ia);
738 continue;
739 }
740 pa = &mesh->edge[ia];
741 pa->tag |= MG_GEO;
742 }
743 }
744 /* get required edges */
745 if ( nedreq ) {
746 rewind(inm);
747 fseek(inm,posnedreq,SEEK_SET);
748 for (k=1; k<=nedreq; k++) {
749 if(!bin) {
750 MMG_FSCANF(inm,"%" MMG5_PRId "",&ia);
751 }
752 else {
753 MMG_FREAD(&ia,MMG5_SW,1,inm);
754 if(iswp) ia=MMG5_swapbin(ia);
755 }
756 if(ia>na) {
757 fprintf(stderr,"\n ## Warning: %s: required Edges number %8" MMG5_PRId "/%8" MMG5_PRId ""
758 " ignored.\n",__func__,ia,na);
759 continue;
760 }
761 pa = &mesh->edge[ia];
762 pa->tag |= MG_REQ;
763 }
764 }
765 /* get parallel edges */
766 if ( nedpar ) {
767 rewind(inm);
768 fseek(inm,posnedpar,SEEK_SET);
769 for (k=1; k<=nedpar; k++) {
770 if(!bin) {
771 MMG_FSCANF(inm,"%" MMG5_PRId "",&ia);
772 }
773 else {
774 MMG_FREAD(&ia,MMG5_SW,1,inm);
775 if(iswp) ia=MMG5_swapbin(ia);
776 }
777 if(ia>na) {
778 fprintf(stderr,"\n ## Warning: %s: parallel Edges number %8" MMG5_PRId "/%8" MMG5_PRId ""
779 " ignored.\n",__func__,ia,na);
780 continue;
781 }
782 pa = &mesh->edge[ia];
783 pa->tag |= MG_PARBDY;
784 }
785 }
786 }
787
788 /* read mesh tetrahedra */
789 rewind(inm);
790 fseek(inm,posne,SEEK_SET);
791 mesh->xt = 0;
792 for (k=1; k<=mesh->ne; k++) {
793 pt = &mesh->tetra[k];
794 if (!bin) {
795 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&pt->v[0],&pt->v[1],&pt->v[2],&pt->v[3],&ref);
796 }
797 else {
798 for (i=0 ; i<4 ; i++) {
799 MMG_FREAD(&pt->v[i],MMG5_SW,1,inm);
800 if(iswp) pt->v[i]=MMG5_swapbin(pt->v[i]);
801 }
802 MMG_FREAD(&ref,MMG5_SW,1,inm);
803 if(iswp) ref=MMG5_swapbin(ref);
804 }
805 if(ref < 0) {
806 nref++;
807 }
808 pt->ref = MMG5_abs(ref);//0;//ref ;
809 for (i=0; i<4; i++) {
810 ppt = &mesh->point[pt->v[i]];
811 ppt->tag &= ~MG_NUL;
812 }
813
814 /* Possibly switch 2 vertex numbers so that each tet is positively oriented */
815 if ( MMG5_orvol(mesh->point,pt->v) < 0.0 ) {
816 /* mesh->xt temporary used to count reoriented tetra */
817 mesh->xt++;
818 aux = pt->v[2];
819 pt->v[2] = pt->v[3];
820 pt->v[3] = aux;
821 }
822 }
823 if(mesh->xt) {
824 fprintf(stderr,"\n $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n");
825 fprintf(stderr," BAD ORIENTATION : vol < 0 -- %8" MMG5_PRId " tetra reoriented\n",mesh->xt);
826 fprintf(stderr," $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n\n");
827 }
828 mesh->xt = 0;
829 /* get required tetrahedra */
830 if(nereq) {
831 rewind(inm);
832 fseek(inm,posnereq,SEEK_SET);
833 for (k=1; k<=nereq; k++) {
834 if(!bin) {
835 MMG_FSCANF(inm,"%d",&i);
836 }
837 else {
838 MMG_FREAD(&i,MMG5_SW,1,inm);
839 if(iswp) i=MMG5_swapbin(i);
840 }
841 if(i>mesh->ne) {
842 fprintf(stderr,"\n ## Warning: %s: required Tetra number %8d"
843 " ignored.\n",__func__,i);
844 continue;
845 }
846 pt = &mesh->tetra[i];
847 pt->tag |= MG_REQ;
848 }
849 }
850 /* get parallel tetrahedra */
851 if(nepar) {
852 rewind(inm);
853 fseek(inm,posnepar,SEEK_SET);
854 for (k=1; k<=nepar; k++) {
855 if(!bin) {
856 MMG_FSCANF(inm,"%d",&i);
857 }
858 else {
859 MMG_FREAD(&i,MMG5_SW,1,inm);
860 if(iswp) i=MMG5_swapbin(i);
861 }
862 if(i>mesh->ne) {
863 fprintf(stderr,"\n ## Warning: %s: parallel Tetra number %8d"
864 " ignored.\n",__func__,i);
865 continue;
866 }
867 pt = &mesh->tetra[i];
868 pt->tag |= MG_PARBDY;
869 }
870 }
871 /* read mesh prisms */
872 rewind(inm);
873 fseek(inm,posnprism,SEEK_SET);
874 for (k=1; k<=mesh->nprism; k++) {
875 pp = &mesh->prism[k];
876 if (!bin) {
877 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&pp->v[0],&pp->v[1],&pp->v[2],
878 &pp->v[3],&pp->v[4],&pp->v[5],&ref);
879 }
880 else {
881 for (i=0 ; i<6 ; i++) {
882 MMG_FREAD(&pp->v[i],MMG5_SW,1,inm);
883 if(iswp) pp->v[i]=MMG5_swapbin(pp->v[i]);
884 }
885 MMG_FREAD(&ref,MMG5_SW,1,inm);
886 if(iswp) ref=MMG5_swapbin(ref);
887 }
888 pp->ref = ref;
889 if ( pp-> ref < 0 ) {
890 pp->ref = -pp->ref;
891 ++nref;
892 }
893 for (i=0; i<6; i++) {
894 ppt = &mesh->point[pp->v[i]];
895 ppt->tag &= ~MG_NUL;
896 }
897 }
898
899 if ( nref ) {
900 fprintf(stdout,"\n $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n");
901 fprintf(stdout," WARNING : %" MMG5_PRId " entities with unexpected refs (ref< 0).\n",nref);
902 fprintf(stdout," We take their absolute values.\n");
903 fprintf(stdout," $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n\n");
904 }
905
906 if ( !mesh->info.iso ) {
907 /* read geometric entities */
908 if ( mesh->nc1 && !ng ) {
909 fprintf(stderr,"\n ## Warning: %s: your mesh don't contains Normals but contains"
910 " NormalAtVertices. The NormalAtVertices are deleted. \n",__func__);
911 mesh->nc1 = 0;
912 }
913
914 if ( ng > 0 ) {
915 MMG5_SAFE_CALLOC(norm,3*ng+1,double,return -1);
916
917 rewind(inm);
918 fseek(inm,posnormal,SEEK_SET);
919 for (k=1; k<=ng; k++) {
920 n = &norm[3*(k-1)+1];
921 if ( mesh->ver == 1 ) {
922 if (!bin) {
923 for (i=0 ; i<3 ; i++) {
924 MMG_FSCANF(inm,"%f",&fc);
925 n[i] = (double) fc;
926 }
927 } else {
928 for (i=0 ; i<3 ; i++) {
929 MMG_FREAD(&fc,MMG5_SW,1,inm);
930 if(iswp) fc=MMG5_swapf(fc);
931 n[i] = (double) fc;
932 }
933 }
934 }
935 else {
936 if (!bin) {
937 MMG_FSCANF(inm,"%lf %lf %lf",&n[0],&n[1],&n[2]);
938 }
939 else {
940 for (i=0 ; i<3 ; i++) {
941 MMG_FREAD(&n[i],MMG5_SD,1,inm);
942 if(iswp) n[i]=MMG5_swapd(n[i]);
943 }
944 }
945 }
946 dd = n[0]*n[0] + n[1]*n[1] + n[2]*n[2];
947 if ( dd > MMG5_EPSD2 ) {
948 dd = 1.0 / sqrt(dd);
949 n[0] *= dd;
950 n[1] *= dd;
951 n[2] *= dd;
952 }
953 }
954
955 rewind(inm);
956 fseek(inm,posnc1,SEEK_SET);
957
958 for (k=1; k<=mesh->nc1; k++) {
959 if (!bin) {
960 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId "",&ip,&idn);
961 }
962 else {
963 MMG_FREAD(&ip,MMG5_SW,1,inm);
964 if(iswp) ip=MMG5_swapbin(ip);
965 MMG_FREAD(&idn,MMG5_SW,1,inm);
966 if(iswp) idn=MMG5_swapbin(idn);
967 }
968 if ( idn > 0 && ip < mesh->np+1 ) {
969 if ( (mesh->info.iso ) && mesh->point[ip].xp == -1 ) {
970 /* Do not store the normals at mesh->info.isoref points (ls mode) */
971 continue;
972 }
973 memcpy(&mesh->point[ip].n,&norm[3*(idn-1)+1],3*sizeof(double));
974 }
975 }
976 MMG5_SAFE_FREE(norm);
977 }
978
979 /* Delete the mark added iso mode */
980 if ( (mesh->info.iso ) && mesh->nc1 ) {
981 for (k=1; k<=mesh->np; k++) {
982 mesh->point[k].xp = 0;
983 }
984 }
985 }
986
987 /* stats */
988 if ( abs(mesh->info.imprim) > 3 ) {
989 fprintf(stdout," NUMBER OF VERTICES %8" MMG5_PRId "\n",mesh->np);
990 fprintf(stdout," NUMBER OF TETRAHEDRA %8" MMG5_PRId "\n",mesh->ne);
991 if ( mesh->nprism )
992 fprintf(stdout," NUMBER OF PRISMS %8" MMG5_PRId "\n",mesh->nprism);
993
994 if ( mesh->na ) {
995 fprintf(stdout," NUMBER OF EDGES %8" MMG5_PRId "\n",mesh->na);
996 if ( nr )
997 fprintf(stdout," NUMBER OF RIDGES %8" MMG5_PRId "\n",nr);
998 }
999 if ( mesh->nt )
1000 fprintf(stdout," NUMBER OF TRIANGLES %8" MMG5_PRId "\n",mesh->nt);
1001 if ( mesh->nquad )
1002 fprintf(stdout," NUMBER OF QUADRILATERALS %8" MMG5_PRId "\n",mesh->nquad);
1003
1004
1005 if ( npreq || nedreq || ntreq || nereq || nqreq ) {
1006 fprintf(stdout," NUMBER OF REQUIRED ENTITIES: \n");
1007 if ( npreq )
1008 fprintf(stdout," VERTICES %8" MMG5_PRId " \n",npreq);
1009 if ( nedreq )
1010 fprintf(stdout," EDGES %8" MMG5_PRId " \n",nedreq);
1011 if ( ntreq )
1012 fprintf(stdout," TRIANGLES %8" MMG5_PRId " \n",ntreq);
1013 if ( nqreq )
1014 fprintf(stdout," QUADRILATERALS %8" MMG5_PRId " \n",nqreq);
1015 if ( nereq )
1016 fprintf(stdout," TETRAHEDRA %8" MMG5_PRId " \n",nereq);
1017 }
1018 if(ncor)
1019 fprintf(stdout," NUMBER OF CORNERS %8" MMG5_PRId " \n",ncor);
1020
1021 if ( nppar || nedpar || ntpar || nepar || nqpar ) {
1022 fprintf(stdout," NUMBER OF PARALLEL ENTITIES: \n");
1023 if ( nppar )
1024 fprintf(stdout," VERTICES %8" MMG5_PRId " \n",nppar);
1025 if ( nedpar )
1026 fprintf(stdout," EDGES %8" MMG5_PRId " \n",nedpar);
1027 if ( ntpar )
1028 fprintf(stdout," TRIANGLES %8" MMG5_PRId " \n",ntpar);
1029 if ( nqpar )
1030 fprintf(stdout," QUADRILATERALS %8" MMG5_PRId " \n",nqpar);
1031 if ( nepar )
1032 fprintf(stdout," TETRAHEDRA %8" MMG5_PRId " \n",nepar);
1033 }
1034 }
1035
1036 return 1;
1037}
1038
1050 FILE* inm;
1051 int bin,ier;
1052
1053 ier = MMG3D_openMesh(mesh->info.imprim,filename,&inm,&bin,"rb","rb");
1054 if( ier < 1 ) return ier;
1055 ier = MMG3D_loadMesh_opened(mesh,inm,bin);
1056 if( ier < 1 ) return ier;
1057
1058 fclose(inm);
1059 return 1;
1060}
1061
1062
1064 FILE* inm;
1065 int ier;
1066 long posNodes,posElts,*posNodeData;
1067 int bin,iswp;
1068 MMG5_int nelts;
1069 int nsols;
1070
1071 mesh->dim = 3;
1072
1074 &posNodes,&posElts,&posNodeData,
1075 &bin,&iswp,&nelts,&nsols);
1076 if ( ier < 1 ) return (ier);
1077
1078 if ( nsols>1 ) {
1079 fprintf(stderr,"Error: SEVERAL SOLUTIONS FOUND (%d)\n",nsols);
1080 fclose(inm);
1081 MMG5_SAFE_FREE(posNodeData);
1082 return -1;
1083 }
1084
1085 if ( !MMG3D_zaldy(mesh) ) {
1086 fclose(inm);
1087 MMG5_SAFE_FREE(posNodeData);
1088 return -1;
1089 }
1090
1091 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt || mesh->nemax < mesh->ne) {
1092 fclose(inm);
1093 MMG5_SAFE_FREE(posNodeData);
1094 return -1;
1095 }
1096
1097 if ( !mesh->ne ) {
1098 fprintf(stderr," ** MISSING DATA.\n");
1099 fprintf(stderr," Check that your mesh contains tetrahedra.\n");
1100 fprintf(stderr," Exit program.\n");
1101 fclose(inm);
1102 MMG5_SAFE_FREE(posNodeData);
1103 return -1;
1104 }
1105
1107 posNodes,posElts,posNodeData,
1108 bin,iswp,nelts,nsols);
1109 MMG5_SAFE_FREE(posNodeData);
1110 if ( ier < 1 ) {
1111 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
1112 return ier;
1113 }
1114
1115 /* Check the metric type */
1116 if ( sol ) {
1117 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,inm);
1118 }
1119
1120 return ier;
1121}
1122
1123/* This is an API function. It is documented in the header.
1124 *
1125 * This function reads a 3D mesh and 0 or 1 data fields in MSH file format (.msh
1126 * extension). We read only low-order points, edges, triangles, quadrangles,
1127 * tetrahedra and prisms.
1128 *
1129 * It returns 0 if the file is not found, -1 if it fails for another reason (mem
1130 * lack, file format...), 1 if success.
1131 */
1133 FILE* inm;
1134 int ier;
1135 long posNodes,posElts,*posNodeData;
1136 int bin,iswp;
1137 MMG5_int nelts;
1138 int nsols;
1139
1140 mesh->dim = 3;
1141
1143 &posNodes,&posElts,&posNodeData,
1144 &bin,&iswp,&nelts,&nsols);
1145 if ( ier < 1 ) return (ier);
1146
1147 mesh->nsols = nsols;
1148 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
1149
1150 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
1151 printf(" Exit program.\n"); fclose(inm);
1152 MMG5_SAFE_FREE(posNodeData);
1153 return -1);
1154 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
1155
1156 if ( !MMG3D_zaldy(mesh) ) {
1157 fclose(inm);
1158 MMG5_SAFE_FREE(posNodeData);
1159 return -1;
1160 }
1161
1162 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt || mesh->nemax < mesh->ne) {
1163 fclose(inm);
1164 MMG5_SAFE_FREE(posNodeData);
1165 return -1;
1166 }
1167
1168 if ( !mesh->ne ) {
1169 fprintf(stderr," ** MISSING DATA.\n");
1170 fprintf(stderr," Check that your mesh contains tetrahedra.\n");
1171 fprintf(stderr," Exit program.\n");
1172 fclose(inm);
1173 MMG5_SAFE_FREE(posNodeData);
1174 return -1;
1175 }
1176
1178 posNodes,posElts,posNodeData,
1179 bin,iswp,nelts,nsols);
1180
1181 if ( ier < 1 ) {
1182 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
1183 }
1184
1185 MMG5_SAFE_FREE(posNodeData);
1186
1187 return ier;
1188}
1189
1191 int ier=0;
1192 const char *filenameptr,*solnameptr;
1193 char *tmp,*soltmp;
1194
1195 if ( filename && strlen(filename) ) {
1196 filenameptr = filename;
1197 solnameptr = filename;
1198 }
1199 else if (mesh->namein && strlen(mesh->namein) ) {
1200 filenameptr = mesh->namein;
1201 if ( sol && strlen(sol->namein) ) {
1202 solnameptr = sol->namein;
1203 }
1204 else {
1205 solnameptr = mesh->namein;
1206 }
1207 }
1208 else {
1209 fprintf(stderr," ## Error: %s: please provide input file name"
1210 " (either in the mesh structure or as function argument).\n",
1211 __func__);
1212 return -1;
1213 }
1214
1215 MMG5_SAFE_MALLOC(tmp,strlen(filenameptr)+1,char,return -1);
1216 strcpy(tmp,filenameptr);
1217
1218 /* read mesh/sol files */
1219 char *ptr = MMG5_Get_filenameExt(tmp);
1220 int fmt = MMG5_Get_format(ptr,MMG5_FMT_MeditASCII);
1221
1222 switch ( fmt ) {
1223
1224 case ( MMG5_FMT_GmshASCII ): case ( MMG5_FMT_GmshBinary ):
1226 break;
1227
1228 case ( MMG5_FMT_VtkVtu ):
1230 break;
1231
1232 case ( MMG5_FMT_VtkVtk ):
1234 break;
1235
1236 case ( MMG5_FMT_MeditASCII ): case ( MMG5_FMT_MeditBinary ):
1238 if ( ier < 1 ) { break; }
1239
1240 /* Optional metric */
1241 if ( sol ) {
1242 MMG5_SAFE_MALLOC(soltmp,strlen(solnameptr)+1,char,return -1);
1243 strcpy(soltmp,solnameptr);
1244
1245 if ( MMG3D_loadSol(mesh,sol,soltmp) == -1) {
1246 fprintf(stderr,"\n ## ERROR: WRONG DATA TYPE OR WRONG SOLUTION NUMBER.\n");
1247 ier = 0;
1248 }
1249 MMG5_SAFE_FREE(soltmp);
1250 }
1251 break;
1252
1253 default:
1254 fprintf(stderr," ** I/O AT FORMAT %s NOT IMPLEMENTED.\n",MMG5_Get_formatName(fmt) );
1255 ier= -1;
1256 }
1257
1259
1260 return ier;
1261}
1262
1274 FILE* inm;
1275 MMG5_pPoint ppt;
1276 MMG5_pTetra pt;
1277 MMG5_pPrism pp;
1278 MMG5_pTria ptt;
1279 MMG5_pQuad pq;
1280 MMG5_xPoint *pxp;
1281 MMG5_int k,na,nc,np,ne,nn,nr,nre,npar,nedreq,nedpar,ntreq,ntpar,nt,nereq,nepar;
1282 MMG5_int npr,nq,nqreq,nqpar,bpos;
1283 int bin,binch;
1284 char chaine[MMG5_FILESTR_LGTH];
1285 static int8_t parWarn = 0;
1286
1287 mesh->ver = 2;
1288
1289 if ( !MMG3D_openMesh(mesh->info.imprim,filename,&inm,&bin,"w","wb") ) {
1290 return 0;
1291 }
1292
1293 /* file header */
1294 binch=0; bpos=10;
1295 if(!bin) {
1296 strcpy(&chaine[0],"MeshVersionFormatted 2\n");
1297 fprintf(inm,"%s",chaine);
1298 strcpy(&chaine[0],"\n\nDimension 3\n");
1299 fprintf(inm,"%s ",chaine);
1300 } else {
1301 binch = 1; //MeshVersionFormatted
1302 fwrite(&binch,MMG5_SW,1,inm);
1303 binch = 2; //version
1304 fwrite(&binch,MMG5_SW,1,inm);
1305 binch = 3; //Dimension
1306 fwrite(&binch,MMG5_SW,1,inm);
1307 bpos = 20; //Pos
1308 fwrite(&bpos,MMG5_SW,1,inm);
1309 binch = 3;
1310 fwrite(&binch,MMG5_SW,1,inm);
1311
1312 }
1313 /* vertices */
1314 np = nc = na = nr = nre = npar = 0;
1315
1316 if ( !mesh->point ) {
1317 fprintf(stderr, "\n ## Error: %s: points array not allocated.\n",
1318 __func__);
1319 fclose(inm);
1320 return 0;
1321 }
1322
1323
1324 for (k=1; k<=mesh->np; k++) {
1325 ppt = &mesh->point[k];
1326 if ( MG_VOK(ppt) ) {
1327 ppt->tmp = ++np;
1328 ppt->flag = 0;
1329 if ( ppt->tag & MG_CRN ) nc++;
1330 if ( ppt->tag & MG_REQ ) nre++;
1331 if ( ppt->tag & MG_PARBDY ) npar++;
1332 }
1333 }
1334
1335 if(!bin) {
1336 strcpy(&chaine[0],"\n\nVertices\n");
1337 fprintf(inm,"%s",chaine);
1338 fprintf(inm,"%" MMG5_PRId "\n",np);
1339 } else {
1340 binch = 4; //Vertices
1341 fwrite(&binch,MMG5_SW,1,inm);
1342 bpos += (3+(1+3*mesh->ver)*np)*MMG5_SW; //NullPos
1343 fwrite(&bpos,MMG5_SW,1,inm);
1344 fwrite(&np,MMG5_SW,1,inm);
1345 }
1346 for (k=1; k<=mesh->np; k++) {
1347 ppt = &mesh->point[k];
1348 if ( MG_VOK(ppt) ) {
1349 if(!bin) {
1350 fprintf(inm,"%.15lg %.15lg %.15lg %" MMG5_PRId "\n",
1351 ppt->c[0],ppt->c[1],ppt->c[2],MMG5_abs(ppt->ref));
1352 } else {
1353 fwrite((unsigned char*)&ppt->c[0],MMG5_SD,1,inm);
1354 fwrite((unsigned char*)&ppt->c[1],MMG5_SD,1,inm);
1355 fwrite((unsigned char*)&ppt->c[2],MMG5_SD,1,inm);
1356 ppt->ref = MMG5_abs(ppt->ref);
1357 fwrite((unsigned char*)&ppt->ref,MMG5_SW,1,inm);
1358 }
1359 }
1360 }
1361
1362 /* corners+required */
1363 if ( nc ) {
1364 if(!bin) {
1365 strcpy(&chaine[0],"\n\nCorners\n");
1366 fprintf(inm,"%s",chaine);
1367 fprintf(inm,"%" MMG5_PRId "\n",nc);
1368 } else {
1369 binch = 13; //
1370 fwrite(&binch,MMG5_SW,1,inm);
1371 bpos += (3+nc)*MMG5_SW; //NullPos
1372 fwrite(&bpos,MMG5_SW,1,inm);
1373 fwrite(&nc,MMG5_SW,1,inm);
1374 }
1375
1376 for (k=1; k<=mesh->np; k++) {
1377 ppt = &mesh->point[k];
1378 if ( MG_VOK(ppt) && ppt->tag & MG_CRN ) {
1379 if(!bin) {
1380 fprintf(inm,"%" MMG5_PRId "\n",ppt->tmp);
1381 } else {
1382 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1383 }
1384 }
1385 }
1386 }
1387 if ( nre ) {
1388 if(!bin) {
1389 strcpy(&chaine[0],"\n\nRequiredVertices\n");
1390 fprintf(inm,"%s",chaine);
1391 fprintf(inm,"%" MMG5_PRId "\n",nre);
1392 } else {
1393 binch = 15; //
1394 fwrite(&binch,MMG5_SW,1,inm);
1395 bpos += (3+nre)*MMG5_SW; //NullPos
1396 fwrite(&bpos,MMG5_SW,1,inm);
1397 fwrite(&nre,MMG5_SW,1,inm);
1398 }
1399 for (k=1; k<=mesh->np; k++) {
1400 ppt = &mesh->point[k];
1401 if ( MG_VOK(ppt) && ppt->tag & MG_REQ ) {
1402 if(!bin) {
1403 fprintf(inm,"%" MMG5_PRId "\n",ppt->tmp);
1404 } else {
1405 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1406 }
1407 }
1408 }
1409 }
1410 if ( npar ) {
1411 if(!bin) {
1412 strcpy(&chaine[0],"\n\nParallelVertices\n");
1413 fprintf(inm,"%s",chaine);
1414 fprintf(inm,"%" MMG5_PRId "\n",npar);
1415 }
1416 else {
1417 if ( !parWarn ) {
1418 parWarn = 1;
1419 fprintf(stderr, "\n ## Warning: %s: parallel entities can't be"
1420 " saved at binary format. Ignored.\n",
1421 __func__);
1422
1423 }
1424 }
1425 for (k=1; k<=mesh->np; k++) {
1426 ppt = &mesh->point[k];
1427 if ( MG_VOK(ppt) && ppt->tag & MG_PARBDY ) {
1428 if(!bin) {
1429 fprintf(inm,"%" MMG5_PRId "\n",ppt->tmp);
1430 }
1431 }
1432 }
1433 }
1434
1435 /* tetrahedra */
1436 ne = nereq = nepar = 0;
1437 if ( mesh->tetra ) {
1438 for (k=1; k<=mesh->ne; k++) {
1439 pt = &mesh->tetra[k];
1440 if ( !MG_EOK(pt) ) {
1441 continue;
1442 }
1443 ne++;
1444 if ( pt->tag & MG_REQ ){
1445 nereq++;
1446 }
1447 if ( pt->tag & MG_PARBDY ){
1448 nepar++;
1449 }
1450 }
1451 }
1452 else {
1453 fprintf(stderr, "\n ## Warning: %s: tetra array not allocated.\n",
1454 __func__);
1455 }
1456
1457 if(!bin) {
1458 strcpy(&chaine[0],"\n\nTetrahedra\n");
1459 fprintf(inm,"%s",chaine);
1460 fprintf(inm,"%" MMG5_PRId "\n",ne);
1461 } else {
1462 binch = 8; //Tetra
1463 fwrite(&binch,MMG5_SW,1,inm);
1464 bpos += (3 + 5*ne)*MMG5_SW;//Pos
1465 fwrite(&bpos,MMG5_SW,1,inm);
1466 fwrite((unsigned char*)&ne,MMG5_SW,1,inm);
1467 }
1468 for (k=1; k<=mesh->ne; k++) {
1469 pt = &mesh->tetra[k];
1470
1471 if ( !MG_EOK(pt) ) { continue; }
1472
1473 /* Tag the tetra vertices to detect points belonging to prisms only (because
1474 * we don't know the normals/tangents at this points, thus we don't want to
1475 * save it). */
1476 mesh->point[pt->v[0]].flag = 1;
1477 mesh->point[pt->v[1]].flag = 1;
1478 mesh->point[pt->v[2]].flag = 1;
1479 mesh->point[pt->v[3]].flag = 1;
1480
1481 if(!bin) {
1482 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",mesh->point[pt->v[0]].tmp,mesh->point[pt->v[1]].tmp
1483 ,mesh->point[pt->v[2]].tmp,mesh->point[pt->v[3]].tmp,pt->ref);
1484 } else {
1485 fwrite(&mesh->point[pt->v[0]].tmp,MMG5_SW,1,inm);
1486 fwrite(&mesh->point[pt->v[1]].tmp,MMG5_SW,1,inm);
1487 fwrite(&mesh->point[pt->v[2]].tmp,MMG5_SW,1,inm);
1488 fwrite(&mesh->point[pt->v[3]].tmp,MMG5_SW,1,inm);
1489 fwrite(&pt->ref,MMG5_SW,1,inm);
1490 }
1491 }
1492
1493 if ( nereq ) {
1494 if(!bin) {
1495 strcpy(&chaine[0],"\n\nRequiredTetrahedra\n");
1496 fprintf(inm,"%s",chaine);
1497 fprintf(inm,"%" MMG5_PRId "\n",nereq);
1498 } else {
1499 binch = 12; //RequiredTetra
1500 fwrite(&binch,MMG5_SW,1,inm);
1501 bpos += (3 + nereq)*MMG5_SW;//Pos
1502 fwrite(&bpos,MMG5_SW,1,inm);
1503 fwrite(&nereq,MMG5_SW,1,inm);
1504 }
1505 ne = 0;
1506 for (k=1; k<=mesh->ne; k++) {
1507 pt = &mesh->tetra[k];
1508 if ( !MG_EOK(pt) ) continue;
1509 ne++;
1510 if ( pt->tag & MG_REQ ) {
1511 if(!bin) {
1512 fprintf(inm,"%" MMG5_PRId "\n",ne);
1513 } else {
1514 fwrite(&ne,MMG5_SW,1,inm);
1515 }
1516 }
1517 }
1518 }
1519
1520 if ( nepar ) {
1521 if(!bin) {
1522 strcpy(&chaine[0],"\n\nParallelTetrahedra\n");
1523 fprintf(inm,"%s",chaine);
1524 fprintf(inm,"%" MMG5_PRId "\n",nepar);
1525 } else {
1526 if ( !parWarn ) {
1527 parWarn = 1;
1528 fprintf(stderr, "\n ## Warning: %s: parallel entities can't be"
1529 " saved at binary format. Ignored.\n",
1530 __func__);
1531
1532 }
1533 }
1534 ne = 0;
1535 for (k=1; k<=mesh->ne; k++) {
1536 pt = &mesh->tetra[k];
1537 if ( !MG_EOK(pt) ) continue;
1538 ne++;
1539 if ( pt->tag & MG_PARBDY ) {
1540 if(!bin) {
1541 fprintf(inm,"%" MMG5_PRId "\n",ne);
1542 }
1543 }
1544 }
1545 }
1546
1547 /* prisms */
1548 npr = 0;
1549 for (k=1; k<=mesh->nprism; k++) {
1550 pp = &mesh->prism[k];
1551 if ( !MG_EOK(pp) ) continue;
1552 npr++;
1553 }
1554
1555 if ( npr ) {
1556 if(!bin) {
1557 strcpy(&chaine[0],"\n\nPrisms\n");
1558 fprintf(inm,"%s",chaine);
1559 fprintf(inm,"%" MMG5_PRId "\n",npr);
1560 } else {
1561 binch = 9; //Prism
1562 fwrite(&binch,MMG5_SW,1,inm);
1563 bpos += (3 + 7*npr)*MMG5_SW;//Pos
1564 fwrite(&bpos,MMG5_SW,1,inm);
1565 fwrite((unsigned char*)&npr,MMG5_SW,1,inm);
1566 }
1567 for (k=1; k<=mesh->nprism; k++) {
1568 pp = &mesh->prism[k];
1569 if ( !MG_EOK(pp) ) continue;
1570
1571 if(!bin) {
1572 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n"
1573 ,mesh->point[pp->v[0]].tmp,mesh->point[pp->v[1]].tmp
1574 ,mesh->point[pp->v[2]].tmp,mesh->point[pp->v[3]].tmp
1575 ,mesh->point[pp->v[4]].tmp,mesh->point[pp->v[5]].tmp,pp->ref);
1576 } else {
1577 fwrite(&mesh->point[pp->v[0]].tmp,MMG5_SW,1,inm);
1578 fwrite(&mesh->point[pp->v[1]].tmp,MMG5_SW,1,inm);
1579 fwrite(&mesh->point[pp->v[2]].tmp,MMG5_SW,1,inm);
1580 fwrite(&mesh->point[pp->v[3]].tmp,MMG5_SW,1,inm);
1581 fwrite(&mesh->point[pp->v[4]].tmp,MMG5_SW,1,inm);
1582 fwrite(&mesh->point[pp->v[5]].tmp,MMG5_SW,1,inm);
1583 fwrite(&pp->ref,MMG5_SW,1,inm);
1584 }
1585 }
1586 }
1587
1588
1589 nn = nt = 0;
1590 if ( mesh->xp && mesh->xpoint ) {
1591 /* Count tangents and normals */
1592 for (k=1; k<=mesh->np; k++) {
1593 ppt = &mesh->point[k];
1594 if ( !MG_VOK(ppt) || (!ppt->flag) || MG_SIN(ppt->tag) ) {
1595 /* We compute normals at some required points but not all of them in analysis
1596 * thus it don't works to save it at the end of the process */
1597 continue;
1598 }
1599 else if ( (ppt->tag & MG_BDY)
1600 && (!(ppt->tag & MG_GEO) || (ppt->tag & MG_NOM)) ) {
1601 /* Indeed non-manifold point along connected mesh have a unique normal if
1602 * pxp->nnor==0 but saving it gives a weird visualization and keeping the
1603 * normal info is useless. */
1604 nn++;
1605 }
1606 if ( MG_EDG_OR_NOM(ppt->tag) ) {
1607 nt++;
1608 }
1609 }
1610
1611 /* write normals */
1612 if(!bin) {
1613 strcpy(&chaine[0],"\n\nNormals\n");
1614 fprintf(inm,"%s",chaine);
1615 fprintf(inm,"%" MMG5_PRId "\n",nn);
1616 } else { binch = 60; //normals
1617 fwrite(&binch,MMG5_SW,1,inm);
1618 bpos += (3+(3*mesh->ver)*nn)*MMG5_SW; //Pos
1619 fwrite(&bpos,MMG5_SW,1,inm);
1620 fwrite(&nn,MMG5_SW,1,inm);
1621 }
1622
1623 for (k=1; k<=mesh->np; k++) {
1624 ppt = &mesh->point[k];
1625 if ( !MG_VOK(ppt) || (!ppt->flag) || MG_SIN(ppt->tag) ) {
1626 continue;
1627 }
1628 else if ( (ppt->tag & MG_BDY)
1629 && (!(ppt->tag & MG_GEO) || (ppt->tag & MG_NOM)) ) {
1630 pxp = &mesh->xpoint[ppt->xp];
1631 if(!bin) {
1632 fprintf(inm,"%.15lg %.15lg %.15lg \n",pxp->n1[0],pxp->n1[1],pxp->n1[2]);
1633 } else {
1634 fwrite((unsigned char*)&pxp->n1[0],MMG5_SD,1,inm);
1635 fwrite((unsigned char*)&pxp->n1[1],MMG5_SD,1,inm);
1636 fwrite((unsigned char*)&pxp->n1[2],MMG5_SD,1,inm);
1637 }
1638 }
1639 }
1640
1641 if(!bin) {
1642 strcpy(&chaine[0],"\n\nNormalAtVertices\n");
1643 fprintf(inm,"%s",chaine);
1644 fprintf(inm,"%" MMG5_PRId "\n",nn);
1645 } else {
1646 binch = 20; //normalatvertices
1647 fwrite(&binch,MMG5_SW,1,inm);
1648 bpos += (3 + 2*nn)*MMG5_SW;//Pos
1649 fwrite(&bpos,MMG5_SW,1,inm);
1650 fwrite(&nn,MMG5_SW,1,inm);
1651 }
1652 nn = 0;
1653 for (k=1; k<=mesh->np; k++) {
1654 ppt = &mesh->point[k];
1655 if ( !MG_VOK(ppt) || (!ppt->flag) || MG_SIN(ppt->tag) ) {
1656 continue;
1657 }
1658 else if ( (ppt->tag & MG_BDY)
1659 && (!(ppt->tag & MG_GEO) || (ppt->tag & MG_NOM) ) ) {
1660 if(!bin) {
1661 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId "\n",ppt->tmp,++nn);
1662 } else {
1663 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1664 ++nn;
1665 fwrite(&nn,MMG5_SW,1,inm);
1666 }
1667 }
1668 }
1669
1670 if ( nt ) {
1671 /* Write tangents */
1672 if(!bin) {
1673 strcpy(&chaine[0],"\n\nTangents\n");
1674 fprintf(inm,"%s",chaine);
1675 fprintf(inm,"%" MMG5_PRId "\n",nt);
1676 } else {
1677 binch = 59; //tangent
1678 fwrite(&binch,MMG5_SW,1,inm);
1679 bpos += (3+(3*mesh->ver)*nt)*MMG5_SW; //Pos
1680 fwrite(&bpos,MMG5_SW,1,inm);
1681 fwrite(&nt,MMG5_SW,1,inm);
1682 }
1683
1684 for (k=1; k<=mesh->np; k++) {
1685 ppt = &mesh->point[k];
1686 if ( !MG_VOK(ppt) || (!ppt->flag) || MG_SIN(ppt->tag) ) {
1687 continue;
1688 }
1689 else if ( MG_EDG_OR_NOM(ppt->tag) ) {
1690 if(!bin) {
1691 fprintf(inm,"%.15lg %.15lg %.15lg \n",ppt->n[0],ppt->n[1],ppt->n[2]);
1692 } else {
1693 fwrite((unsigned char*)&ppt->n[0],MMG5_SD,1,inm);
1694 fwrite((unsigned char*)&ppt->n[1],MMG5_SD,1,inm);
1695 fwrite((unsigned char*)&ppt->n[2],MMG5_SD,1,inm);
1696 }
1697 }
1698 }
1699
1700
1701 if(!bin) {
1702 strcpy(&chaine[0],"\n\nTangentAtVertices\n");
1703 fprintf(inm,"%s",chaine);
1704 fprintf(inm,"%" MMG5_PRId "\n",nt);
1705 } else {
1706 binch = 61; //tangentatvertices
1707 fwrite(&binch,MMG5_SW,1,inm);
1708 bpos += (3 + 2*nt)*MMG5_SW;//Pos
1709 fwrite(&bpos,MMG5_SW,1,inm);
1710 fwrite(&nt,MMG5_SW,1,inm);
1711 }
1712 nt = 0;
1713 for (k=1; k<=mesh->np; k++) {
1714 ppt = &mesh->point[k];
1715 if ( !MG_VOK(ppt) || (!ppt->flag) || MG_SIN(ppt->tag) ) {
1716 continue;
1717 }
1718 else if ( MG_EDG_OR_NOM(ppt->tag) ) {
1719 if(!bin) {
1720 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId "\n",ppt->tmp,++nt);
1721 } else {
1722 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1723 ++nt;
1724 fwrite(&(nn),MMG5_SW,1,inm);
1725 }
1726 }
1727 }
1728 }
1729 }
1730
1731 /* boundary mesh */
1732 /* tria + required tria */
1733 ntreq = ntpar = 0;
1734
1735 if ( mesh->nt ) {
1736 if(!bin) {
1737 strcpy(&chaine[0],"\n\nTriangles\n");
1738 fprintf(inm,"%s",chaine);
1739 fprintf(inm,"%" MMG5_PRId "\n",mesh->nt);
1740 } else {
1741 binch = 6; //Triangles
1742 fwrite(&binch,MMG5_SW,1,inm);
1743 bpos += (3+4*mesh->nt)*MMG5_SW; //Pos
1744 fwrite(&bpos,MMG5_SW,1,inm);
1745 fwrite(&mesh->nt,MMG5_SW,1,inm);
1746 }
1747 for (k=1; k<=mesh->nt; k++) {
1748 ptt = &mesh->tria[k];
1749 if ( ptt->tag[0] & MG_REQ && ptt->tag[1] & MG_REQ && ptt->tag[2] & MG_REQ ) {
1750 ntreq++;
1751 }
1752 if ( ptt->tag[0] & MG_PARBDY && ptt->tag[1] & MG_PARBDY && ptt->tag[2] & MG_PARBDY ) {
1753 ntpar++;
1754 }
1755 if(!bin) {
1756 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",mesh->point[ptt->v[0]].tmp,mesh->point[ptt->v[1]].tmp
1757 ,mesh->point[ptt->v[2]].tmp,ptt->ref);
1758 } else {
1759 fwrite(&mesh->point[ptt->v[0]].tmp,MMG5_SW,1,inm);
1760 fwrite(&mesh->point[ptt->v[1]].tmp,MMG5_SW,1,inm);
1761 fwrite(&mesh->point[ptt->v[2]].tmp,MMG5_SW,1,inm);
1762 fwrite(&ptt->ref,MMG5_SW,1,inm);
1763 }
1764 }
1765 if ( ntreq ) {
1766 if(!bin) {
1767 strcpy(&chaine[0],"\n\nRequiredTriangles\n");
1768 fprintf(inm,"%s",chaine);
1769 fprintf(inm,"%" MMG5_PRId "\n",ntreq);
1770 } else {
1771 binch = 17; //ReqTriangles
1772 fwrite(&binch,MMG5_SW,1,inm);
1773 bpos += (3+ntreq)*MMG5_SW; //Pos
1774 fwrite(&bpos,MMG5_SW,1,inm);
1775 fwrite(&ntreq,MMG5_SW,1,inm);
1776 }
1777 for (k=0; k<=mesh->nt; k++) {
1778 ptt = &mesh->tria[k];
1779 if ( (ptt->tag[0] & MG_REQ) && (ptt->tag[1] & MG_REQ)
1780 && ptt->tag[2] & MG_REQ ) {
1781 if(!bin) {
1782 fprintf(inm,"%" MMG5_PRId "\n",k);
1783 } else {
1784 fwrite(&k,MMG5_SW,1,inm);
1785 }
1786 }
1787 }
1788 }
1789 if ( ntpar ) {
1790 if(!bin) {
1791 strcpy(&chaine[0],"\n\nParallelTriangles\n");
1792 fprintf(inm,"%s",chaine);
1793 fprintf(inm,"%" MMG5_PRId "\n",ntpar);
1794 } else {
1795 if ( !parWarn ) {
1796 parWarn = 1;
1797 fprintf(stderr, "\n ## Warning: %s: parallel entities can't be"
1798 " saved at binary format. Ignored.\n",
1799 __func__);
1800 }
1801 }
1802 for (k=0; k<=mesh->nt; k++) {
1803 ptt = &mesh->tria[k];
1804 if ( (ptt->tag[0] & MG_PARBDY) && (ptt->tag[1] & MG_PARBDY)
1805 && ptt->tag[2] & MG_PARBDY ) {
1806 if(!bin) {
1807 fprintf(inm,"%" MMG5_PRId "\n",k);
1808 }
1809 }
1810 }
1811 }
1812 }
1813
1814 /* quad + required quad */
1815 nq = nqreq = nqpar = 0;
1816
1817 if ( mesh->nquad ) {
1818
1819 for (k=1; k<=mesh->nquad; k++) {
1820 pq = &mesh->quadra[k];
1821 if ( !MG_EOK(pq) ) {
1822 continue;
1823 }
1824 nq++;
1825 if ( pq->tag[0] & MG_REQ && pq->tag[1] & MG_REQ &&
1826 pq->tag[2] & MG_REQ && pq->tag[3] & MG_REQ ) {
1827 nqreq++;
1828 }
1829 if ( pq->tag[0] & MG_PARBDY && pq->tag[1] & MG_PARBDY &&
1830 pq->tag[2] & MG_PARBDY && pq->tag[3] & MG_PARBDY ) {
1831 nqpar++;
1832 }
1833 }
1834 }
1835
1836 if ( nq ) {
1837 if(!bin) {
1838 strcpy(&chaine[0],"\n\nQuadrilaterals\n");
1839 fprintf(inm,"%s",chaine);
1840 fprintf(inm,"%" MMG5_PRId "\n",nq);
1841 } else {
1842 binch = 7; //Quadrilaterals
1843 fwrite(&binch,MMG5_SW,1,inm);
1844 bpos += (3+5*nq)*MMG5_SW; //Pos
1845 fwrite(&bpos,MMG5_SW,1,inm);
1846 fwrite(&nq,MMG5_SW,1,inm);
1847 }
1848 for (k=1; k<=mesh->nquad; k++) {
1849 pq = &mesh->quadra[k];
1850 if ( !MG_EOK(pq) ) continue;
1851
1852 if(!bin) {
1853 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",mesh->point[pq->v[0]].tmp,
1854 mesh->point[pq->v[1]].tmp,mesh->point[pq->v[2]].tmp,
1855 mesh->point[pq->v[3]].tmp, pq->ref);
1856 } else {
1857 fwrite(&mesh->point[pq->v[0]].tmp,MMG5_SW,1,inm);
1858 fwrite(&mesh->point[pq->v[1]].tmp,MMG5_SW,1,inm);
1859 fwrite(&mesh->point[pq->v[2]].tmp,MMG5_SW,1,inm);
1860 fwrite(&mesh->point[pq->v[3]].tmp,MMG5_SW,1,inm);
1861 fwrite(&pq->ref,MMG5_SW,1,inm);
1862 }
1863 }
1864 if ( nqreq ) {
1865 if(!bin) {
1866 strcpy(&chaine[0],"\n\nRequiredQuadrilaterals\n");
1867 fprintf(inm,"%s",chaine);
1868 fprintf(inm,"%" MMG5_PRId "\n",nqreq);
1869 } else {
1870 binch = 18; //ReqQuad
1871 fwrite(&binch,MMG5_SW,1,inm);
1872 bpos += (3+nqreq)*MMG5_SW; //Pos
1873 fwrite(&bpos,MMG5_SW,1,inm);
1874 fwrite(&nqreq,MMG5_SW,1,inm);
1875 }
1876 for (k=0; k<=mesh->nquad; k++) {
1877 pq = &mesh->quadra[k];
1878 if ( (pq->tag[0] & MG_REQ) && (pq->tag[1] & MG_REQ)
1879 && pq->tag[2] & MG_REQ && pq->tag[3] & MG_REQ ) {
1880 if(!bin) {
1881 fprintf(inm,"%" MMG5_PRId "\n",k);
1882 } else {
1883 fwrite(&k,MMG5_SW,1,inm);
1884 }
1885 }
1886 }
1887 }
1888 if ( nqpar ) {
1889 if(!bin) {
1890 strcpy(&chaine[0],"\n\nParallelQuadrilaterals\n");
1891 fprintf(inm,"%s",chaine);
1892 fprintf(inm,"%" MMG5_PRId "\n",nqpar);
1893 } else {
1894 if ( !parWarn ) {
1895 parWarn = 1;
1896 fprintf(stderr, "\n ## Warning: %s: parallel entities can't be"
1897 " saved at binary format. Ignored.\n",
1898 __func__);
1899 }
1900 }
1901 for (k=0; k<=mesh->nquad; k++) {
1902 pq = &mesh->quadra[k];
1903 if ( (pq->tag[0] & MG_PARBDY) && (pq->tag[1] & MG_PARBDY) &&
1904 pq->tag[2] & MG_PARBDY && pq->tag[3] & MG_PARBDY ) {
1905 if(!bin) {
1906 fprintf(inm,"%" MMG5_PRId "\n",k);
1907 }
1908 }
1909 }
1910 }
1911 }
1912
1913 nr = nedreq = nedpar = 0;
1914 if ( mesh->na ) {
1915 if(!bin) {
1916 strcpy(&chaine[0],"\n\nEdges\n");
1917 fprintf(inm,"%s",chaine);
1918 fprintf(inm,"%" MMG5_PRId "\n",mesh->na);
1919 } else {
1920 binch = 5; //Edges
1921 fwrite(&binch,MMG5_SW,1,inm);
1922 bpos += (3 + 3*mesh->na)*MMG5_SW;//Pos
1923 fwrite(&bpos,MMG5_SW,1,inm);
1924 fwrite(&mesh->na,MMG5_SW,1,inm);
1925 }
1926 for (k=1; k<=mesh->na; k++) {
1927 if(!bin) {
1928 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",mesh->point[mesh->edge[k].a].tmp,
1929 mesh->point[mesh->edge[k].b].tmp,mesh->edge[k].ref);
1930 } else {
1931 fwrite(&mesh->point[mesh->edge[k].a].tmp,MMG5_SW,1,inm);
1932 fwrite(&mesh->point[mesh->edge[k].b].tmp,MMG5_SW,1,inm);
1933 fwrite(&mesh->edge[k].ref,MMG5_SW,1,inm);
1934 }
1935 if ( mesh->edge[k].tag & MG_GEO ) nr++;
1936 if ( mesh->edge[k].tag & MG_REQ ) nedreq++;
1937 if ( mesh->edge[k].tag & MG_PARBDY ) nedpar++;
1938 }
1939
1940 if ( nr ) {
1941 if(!bin) {
1942 strcpy(&chaine[0],"\n\nRidges\n");
1943 fprintf(inm,"%s",chaine);
1944 fprintf(inm,"%" MMG5_PRId "\n",nr);
1945 } else {
1946 binch = 14; //Ridges
1947 fwrite(&binch,MMG5_SW,1,inm);
1948 bpos += (3 + nr)*MMG5_SW;//Pos
1949 fwrite(&bpos,MMG5_SW,1,inm);
1950 fwrite(&nr,MMG5_SW,1,inm);
1951 }
1952 na = 0;
1953 for (k=1; k<=mesh->na; k++) {
1954 na++;
1955 if ( mesh->edge[k].tag & MG_GEO ) {
1956 if(!bin) {
1957 fprintf(inm,"%" MMG5_PRId "\n",na);
1958 } else {
1959 fwrite(&na,MMG5_SW,1,inm);
1960 }
1961 }
1962 }
1963 }
1964
1965 if ( nedreq ) {
1966 if(!bin) {
1967 strcpy(&chaine[0],"\n\nRequiredEdges\n");
1968 fprintf(inm,"%s",chaine);
1969 fprintf(inm,"%" MMG5_PRId "\n",nedreq);
1970 } else {
1971 binch = 16; //RequiredEdges
1972 fwrite(&binch,MMG5_SW,1,inm);
1973 bpos += (3 + nedreq)*MMG5_SW;//Pos
1974 fwrite(&bpos,MMG5_SW,1,inm);
1975 fwrite(&nedreq,MMG5_SW,1,inm);
1976 }
1977 na = 0;
1978 for (k=1; k<=mesh->na; k++) {
1979 na++;
1980 if ( mesh->edge[k].tag & MG_REQ ) {
1981 if(!bin) {
1982 fprintf(inm,"%" MMG5_PRId "\n",na);
1983 } else {
1984 fwrite(&na,MMG5_SW,1,inm);
1985 }
1986 }
1987 }
1988 }
1989 if ( nedpar ) {
1990 if(!bin) {
1991 strcpy(&chaine[0],"\n\nParallelEdges\n");
1992 fprintf(inm,"%s",chaine);
1993 fprintf(inm,"%" MMG5_PRId "\n",nedpar);
1994 } else {
1995 if ( !parWarn ) {
1996 parWarn = 1;
1997 fprintf(stderr, "\n ## Warning: %s: parallel entities can't be"
1998 " saved at binary format. Ignored.\n",
1999 __func__);
2000 }
2001 }
2002 na = 0;
2003 for (k=1; k<=mesh->na; k++) {
2004 na++;
2005 if ( mesh->edge[k].tag & MG_PARBDY ) {
2006 if(!bin) {
2007 fprintf(inm,"%" MMG5_PRId "\n",na);
2008 }
2009 }
2010 }
2011 }
2012 }
2013
2014
2015 if ( mesh->info.imprim > 4 ) {
2016 fprintf(stdout," NUMBER OF VERTICES %8" MMG5_PRId " CORNERS %8" MMG5_PRId ""
2017 " REQUIRED %8" MMG5_PRId "\n",np,nc,nre);
2018 fprintf(stdout," NUMBER OF TETRAHEDRA %8" MMG5_PRId " REQUIRED %8" MMG5_PRId "\n",
2019 ne,nereq);
2020 if ( npr )
2021 fprintf(stdout," NUMBER OF PRISMS %8" MMG5_PRId "\n",npr);
2022
2023 if ( na )
2024 fprintf(stdout," NUMBER OF EDGES %8" MMG5_PRId " RIDGES %8" MMG5_PRId ""
2025 " REQUIRED %8" MMG5_PRId "\n",na,nr,nedreq);
2026 if ( mesh->nt )
2027 fprintf(stdout," NUMBER OF TRIANGLES %8" MMG5_PRId " REQUIRED %8" MMG5_PRId "\n",
2028 mesh->nt, ntreq);
2029 if ( nq )
2030 fprintf(stdout," NUMBER OF QUADRILATERALS %8" MMG5_PRId " REQUIRED %8" MMG5_PRId "\n",
2031 nq,nqreq);
2032
2033 if ( npar || nedpar || ntpar || nepar || nqpar ) {
2034 fprintf(stdout," NUMBER OF PARALLEL ENTITIES: \n");
2035 if ( npar )
2036 fprintf(stdout," VERTICES %8" MMG5_PRId " \n",npar);
2037 if ( nedpar )
2038 fprintf(stdout," EDGES %8" MMG5_PRId " \n",nedpar);
2039 if ( ntpar )
2040 fprintf(stdout," TRIANGLES %8" MMG5_PRId " \n",ntpar);
2041 if ( nqpar )
2042 fprintf(stdout," QUADRILATERALS %8" MMG5_PRId " \n",nqpar);
2043 if ( nepar )
2044 fprintf(stdout," TETRAHEDRA %8" MMG5_PRId " \n",nepar);
2045 }
2046 }
2047
2048 /* end of file */
2049 if(!bin) {
2050 strcpy(&chaine[0],"\n\nEnd\n");
2051 fprintf(inm,"%s",chaine);
2052 } else {
2053 binch = 54; //End
2054 fwrite(&binch,MMG5_SW,1,inm);
2055 bpos += 2*MMG5_SW; //bpos + End key
2056 fwrite(&bpos,MMG5_SW,1,inm);
2057 }
2058 fclose(inm);
2059 return 1;
2060}
2061
2063 int ier=0;
2064 const char *filenameptr,*solnameptr;
2065 char *tmp,*soltmp;
2066
2067 if ( filename && strlen(filename) ) {
2068 filenameptr = filename;
2069 solnameptr = filename;
2070 }
2071 else if (mesh->namein && strlen(mesh->namein) ) {
2072 filenameptr = mesh->namein;
2073 if ( sol && strlen(sol->namein) ) {
2074 solnameptr = sol->namein;
2075 }
2076 else {
2077 solnameptr = mesh->namein;
2078 }
2079 }
2080 else {
2081 fprintf(stderr," ## Error: %s: please provide input file name"
2082 " (either in the mesh structure or as function argument).\n",
2083 __func__);
2084 return 0;
2085 }
2086
2087 MMG5_SAFE_MALLOC(tmp,strlen(filenameptr)+1,char,return 0);
2088 strcpy(tmp,filenameptr);
2089
2090 /* read mesh/sol files */
2091 char *ptr = MMG5_Get_filenameExt(tmp);
2092 int fmt = MMG5_Get_format(ptr,MMG5_FMT_MeditASCII);
2093
2094 int8_t savesolFile = 0;
2095
2096 switch ( fmt ) {
2097 case ( MMG5_FMT_GmshASCII ): case ( MMG5_FMT_GmshBinary ):
2099 break;
2100 case ( MMG5_FMT_VtkVtu ):
2102 break;
2103 case ( MMG5_FMT_VtkVtk ):
2105 break;
2106 case ( MMG5_FMT_Tetgen ):
2108 savesolFile = 1;
2109 break;
2110 default:
2112 savesolFile = 1;
2113 break;
2114 }
2115
2116 if ( ier && savesolFile ) {
2117 /* Medit or tetgen output: save the solution in a .sol file */
2118 if ( sol && sol->np ) {
2119 MMG5_SAFE_MALLOC(soltmp,strlen(solnameptr)+1,char,return 0);
2120 strcpy(soltmp,solnameptr);
2121
2122 if ( MMG3D_saveSol(mesh,sol,soltmp) == -1) {
2123 fprintf(stderr,"\n ## ERROR: WRONG DATA TYPE OR WRONG SOLUTION NUMBER.\n");
2124 ier = 0;
2125 }
2126 MMG5_SAFE_FREE(soltmp);
2127 }
2128 }
2129
2130 return ier;
2131}
2132
2134
2135 return MMG5_saveMshMesh(mesh,&sol,filename,1);
2136}
2137
2139
2141}
2142
2144 FILE *inm;
2145 long posnp;
2146 int iswp,ier,ver,bin,*type,nsols,dim;
2147 MMG5_int k,np;
2148
2150 ier = MMG5_loadSolHeader(filename,3,&inm,&ver,&bin,&iswp,&np,&dim,&nsols,
2151 &type,&posnp,mesh->info.imprim);
2152 if ( ier < 1 ) return ier;
2153
2154 if ( nsols!=1 ) {
2155 fprintf(stderr,"Error: SEVERAL SOLUTIONS FOUND (%d)\n",nsols);
2156 fclose(inm);
2157 MMG5_SAFE_FREE(type);
2158 return -1;
2159 }
2160
2161 if ( mesh->np != np ) {
2162 fprintf(stderr," ** MISMATCHES DATA: THE NUMBER OF VERTICES IN "
2163 "THE MESH (%" MMG5_PRId ") DIFFERS FROM THE NUMBER OF VERTICES IN "
2164 "THE SOLUTION (%" MMG5_PRId ") \n",mesh->np,np);
2165 fclose(inm);
2166 MMG5_SAFE_FREE(type);
2167 return -1;
2168 }
2169
2170 /* #MMG5_loadSolHeader function reads only solution at vertices so we don't
2171 have to check the entites on which the metric applies */
2172 int entities = MMG5_Vertex;
2173 ier = MMG5_chkMetricType(mesh,type,&entities,inm);
2174
2175 if ( ier < 1 ) {
2176 MMG5_SAFE_FREE(type);
2177 return ier;
2178 }
2179
2180 /* Allocate and store the header information for each solution */
2181 if ( !MMG3D_Set_solSize(mesh,met,MMG5_Vertex,mesh->np,type[0]) ) {
2182 fclose(inm);
2183 MMG5_SAFE_FREE(type);
2184 return -1;
2185 }
2186 /* For binary files, we read the verson inside the file */
2187 if ( ver ) met->ver = ver;
2188
2189 MMG5_SAFE_FREE(type);
2190
2191 /* Read mesh solutions */
2192 rewind(inm);
2193 fseek(inm,posnp,SEEK_SET);
2194
2195 if ( met->ver == 1 ) {
2196 /* Single precision */
2197 for (k=1; k<=mesh->np; k++) {
2198 if ( MMG5_readFloatSol3D(met,inm,bin,iswp,k) < 0 ) return -1;
2199 }
2200 }
2201 else {
2202 /* Double precision */
2203 for (k=1; k<=mesh->np; k++) {
2204 if ( MMG5_readDoubleSol3D(met,inm,bin,iswp,k) < 0 ) return -1;
2205 }
2206 }
2207
2208 fclose(inm);
2209
2210 /* stats */
2212
2213 return 1;
2214}
2215
2217 MMG5_pSol psl;
2218 FILE *inm;
2219 long posnp;
2220 int iswp,ier,ver,bin,*type,nsols,dim;
2221 MMG5_int j,k,np;
2222 char data[16];
2223 static char mmgWarn = 0;
2224
2226 ier = MMG5_loadSolHeader(filename,3,&inm,&ver,&bin,&iswp,&np,&dim,&nsols,
2227 &type,&posnp,mesh->info.imprim);
2228 if ( ier < 1 ) return ier;
2229
2230 if ( mesh->np != np ) {
2231 fprintf(stderr," ** MISMATCHES DATA: THE NUMBER OF VERTICES IN "
2232 "THE MESH (%" MMG5_PRId ") DIFFERS FROM THE NUMBER OF VERTICES IN "
2233 "THE SOLUTION (%" MMG5_PRId ") \n",mesh->np,np);
2234 fclose(inm);
2235 MMG5_SAFE_FREE(type);
2236 return -1;
2237 }
2238
2240 mesh->nsols = nsols;
2241
2242 if ( nsols > MMG5_NSOLS_MAX ) {
2243 fprintf(stderr,"\n ## Error: %s: unexpected number of data (%d).\n",
2244 __func__,nsols);
2245 MMG5_SAFE_FREE(type);
2246 fclose(inm);
2247 return -1;
2248 }
2249
2250 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
2251
2252 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
2253 printf(" Exit program.\n"); fclose(inm);
2254 MMG5_SAFE_FREE(type);
2255 return -1);
2256 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
2257
2258 for ( j=0; j<nsols; ++j ) {
2259 psl = *sol + j;
2260
2261 /* Give an arbitrary name to the solution because the Medit format has no
2262 * name field */
2263 sprintf(data,"sol_%" MMG5_PRId "",j);
2264 if ( !MMG3D_Set_inputSolName(mesh,psl,data) ) {
2265 if ( !mmgWarn ) {
2266 mmgWarn = 1;
2267 fprintf(stderr,"\n ## Warning: %s: unable to set solution name for"
2268 " at least 1 solution.\n",__func__);
2269 }
2270 }
2271
2272 /* Allocate and store the header informations for each solution */
2273 if ( !MMG3D_Set_solSize(mesh,psl,MMG5_Vertex,mesh->np,type[j]) ) {
2274 MMG5_SAFE_FREE(type);
2275 fclose(inm);
2276 return -1;
2277 }
2278
2279 /* For binary file, we read the version inside the file */
2280 if ( ver ) psl->ver = ver;
2281 }
2282 MMG5_SAFE_FREE(type);
2283
2284 /* read mesh solutions */
2285 rewind(inm);
2286 fseek(inm,posnp,SEEK_SET);
2287
2288 if ( (*sol)[0].ver == 1 ) {
2289 /* Single precision */
2290 for (k=1; k<=mesh->np; k++) {
2291 for ( j=0; j<nsols; ++j ) {
2292 psl = *sol + j;
2293 if ( MMG5_readFloatSol3D(psl,inm,bin,iswp,k) < 0 ) return -1;
2294 }
2295 }
2296 }
2297 else {
2298 /* Double precision */
2299 for (k=1; k<=mesh->np; k++) {
2300 for ( j=0; j<nsols; ++j ) {
2301 psl = *sol + j;
2302 if ( MMG5_readDoubleSol3D(psl,inm,bin,iswp,k) < 0 ) return -1;
2303 }
2304 }
2305 }
2306 fclose(inm);
2307
2308 /* stats */
2310
2311 return 1;
2312}
2313
2315 FILE* inm;
2316 MMG5_pPoint ppt;
2317 int binch,bin,ier;
2318 MMG5_int k,bpos;
2319
2320 if ( !met->m ) {
2321 fprintf(stderr,"\n ## Warning: %s: no metric data to save.\n",__func__);
2322 return 1;
2323 }
2324
2325 met->ver = 2;
2326
2327 bpos = 0;
2328 ier = MMG5_saveSolHeader( mesh,filename,&inm,met->ver,&bin,&bpos,mesh->np,
2329 met->dim,1,&met->entities,&met->type,&met->size);
2330
2331 if ( ier < 1 ) return ier;
2332
2333 for (k=1; k<=mesh->np; k++) {
2334 ppt = &mesh->point[k];
2335 if ( !MG_VOK(ppt) ) continue;
2336
2337 MMG5_writeDoubleSol3D(mesh,met,inm,bin,k,1);
2338 fprintf(inm,"\n");
2339 }
2340
2341 /* End file */
2342 if(!bin) {
2343 fprintf(inm,"\n\nEnd\n");
2344 } else {
2345 binch = 54; //End
2346 fwrite(&binch,MMG5_SW,1,inm);
2347 }
2348 fclose(inm);
2349 return 1;
2350}
2351
2353 MMG5_pSol psl;
2354 FILE* inm;
2355 MMG5_pPoint ppt;
2356 MMG5_pTetra pt;
2357 int binch,bin,ier,npointSols,ncellSols;
2358 MMG5_int j,k,bpos;
2359 int *type,*entities,*size;
2360
2361 if ( !(*sol)[0].m ) return -1;
2362
2363 (*sol)[0].ver = 2;
2364
2365 MMG5_SAFE_CALLOC(type,mesh->nsols,int,return 0);
2366 MMG5_SAFE_CALLOC(size,mesh->nsols,int,MMG5_SAFE_FREE(type);return 0);
2367 MMG5_SAFE_CALLOC(entities,mesh->nsols,int,
2368 MMG5_SAFE_FREE(type);MMG5_SAFE_FREE(size);return 0);
2369
2370 npointSols = 0;
2371 ncellSols = 0;
2372 for (k=0; k<mesh->nsols; ++k ) {
2373 if ( ((*sol)[k].entities==MMG5_Noentity) || ((*sol)[k].entities==MMG5_Vertex) ) {
2374 ++npointSols;
2375 }
2376 else if ( (*sol)[k].entities == MMG5_Tetrahedron ) {
2377 ++ncellSols;
2378 }
2379 else {
2380 printf("\n ## Warning: %s: unexpected entity type for solution %" MMG5_PRId ": %s."
2381 "\n Ignored.\n",
2382 __func__,k,MMG5_Get_entitiesName((*sol)[k].entities));
2383 }
2384 type[k] = (*sol)[k].type;
2385 size[k] = (*sol)[k].size;
2386 entities[k] = (*sol)[k].entities;
2387 }
2388
2389 bpos = 0;
2390 ier = MMG5_saveSolHeader( mesh,filename,&inm,(*sol)[0].ver,&bin,&bpos,mesh->np,
2391 (*sol)[0].dim,mesh->nsols,entities,type,size);
2392
2393 if ( ier < 1 ) return ier;
2394
2395 for (k=1; k<=mesh->np; k++) {
2396 ppt = &mesh->point[k];
2397 if ( !MG_VOK(ppt) ) continue;
2398
2399 for ( j=0; j<mesh->nsols; ++j ) {
2400 psl = *sol+j;
2401
2402 if ( (psl->entities==MMG5_Noentity) || (psl->entities==MMG5_Vertex) ) {
2403 MMG5_writeDoubleSol3D(mesh,psl,inm,bin,k,0);
2404 }
2405 }
2406 fprintf(inm,"\n");
2407 }
2408
2409 MMG5_saveSolAtTetrahedraHeader( mesh,inm,(*sol)[0].ver,bin,&bpos,mesh->nsols,
2410 ncellSols,entities,type,size );
2411
2412 for (k=1; k<=mesh->ne; k++) {
2413 pt = &mesh->tetra[k];
2414 if ( !MG_EOK(pt) ) continue;
2415
2416 for ( j=0; j<mesh->nsols; ++j ) {
2417 psl = *sol+j;
2418 if ( psl->entities==MMG5_Tetrahedron ) {
2419 MMG5_writeDoubleSol3D(mesh,psl,inm,bin,k,0);
2420 }
2421 }
2422 fprintf(inm,"\n");
2423 }
2424
2425
2426 MMG5_SAFE_FREE(type);
2427 MMG5_SAFE_FREE(size);
2428 MMG5_SAFE_FREE(entities);
2429
2430 /* End file */
2431 if(!bin) {
2432 fprintf(inm,"\n\nEnd\n");
2433 } else {
2434 binch = 54; //End
2435 fwrite(&binch,MMG5_SW,1,inm);
2436 }
2437 fclose(inm);
2438 return 1;
2439}
2440
2441static inline
2443 FILE* inm;
2444 MMG5_pTetra pt;
2445 MMG5_int k,i,ne;
2446 char *ptr,*data;
2447
2448 if ( !mesh->ne ) {
2449 return 1;
2450 }
2451
2452 if ( (!filename) || !(*filename) ) {
2454 }
2455 if ( (!filename) || !(*filename) ) {
2456 printf("\n ## Error: %s: unable to save a file without a valid filename\n.",
2457 __func__);
2458 return 0;
2459 }
2460
2461 /* Name of file */
2462 MMG5_SAFE_CALLOC(data,strlen(filename)+5,char,return 0);
2463 strcpy(data,filename);
2464 ptr = strstr(data,".node");
2465 if ( ptr ) {
2466 *ptr = '\0';
2467 }
2468
2469 /* Add .node ext */
2470 strcat(data,".ele");
2471 if( !(inm = fopen(data,"wb")) ) {
2472 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
2473 MMG5_SAFE_FREE(data);
2474 return 0;
2475 }
2476
2477 fprintf(stdout," %%%% %s OPENED\n",data);
2478 MMG5_SAFE_FREE(data);
2479
2480 ne = 0;
2481 for (k=1; k<=mesh->ne; k++) {
2482 pt = &mesh->tetra[k];
2483 if ( !MG_EOK(pt) ) continue;
2484 ne++;
2485 }
2486
2487 /* Save elt number, node number per elt, 1 bdy marker */
2488 fprintf(inm, "%" MMG5_PRId " %d %d\n\n",ne,mesh->dim+1,1);
2489
2490 ne = 0;
2491 for ( k=1; k<=mesh->ne; ++k ) {
2492 pt = &mesh->tetra[k];
2493 if ( MG_EOK(pt) ) {
2494 /* Save elt idx */
2495 fprintf(inm, "%" MMG5_PRId " ",++ne);
2496
2497 /* Save connectivity */
2498 for ( i=0; i<=mesh->dim; ++i ) {
2499 fprintf(inm, "%" MMG5_PRId " ",mesh->point[pt->v[i]].tmp);
2500 }
2501
2502 /* Save bdy marker */
2503 fprintf(inm, "%" MMG5_PRId "\n",pt->ref);
2504 }
2505 }
2506 fprintf(stdout," NUMBER OF ELEMENT %8" MMG5_PRId "\n",ne);
2507
2508 fclose(inm);
2509
2510 return 1;
2511}
2512
2513static inline
2515 FILE* inm;
2516 MMG5_pTetra pt;
2517 MMG5_int k,i,ne;
2518 MMG5_int idx;
2519 char *ptr,*data;
2520
2521 if ( !mesh->na ) {
2522 return 1;
2523 }
2524
2525 if ( (!filename) || !(*filename) ) {
2527 }
2528 if ( (!filename) || !(*filename) ) {
2529 printf("\n ## Error: %s: unable to save a file without a valid filename\n.",
2530 __func__);
2531 return 0;
2532 }
2533
2534 /* Name of file */
2535 MMG5_SAFE_CALLOC(data,strlen(filename)+7,char,return 0);
2536 strcpy(data,filename);
2537 ptr = strstr(data,".node");
2538 if ( ptr ) {
2539 *ptr = '\0';
2540 }
2541
2542 /* Add .node ext */
2543 strcat(data,".neigh");
2544 if( !(inm = fopen(data,"wb")) ) {
2545 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
2546 MMG5_SAFE_FREE(data);
2547 return 0;
2548 }
2549
2550 fprintf(stdout," %%%% %s OPENED\n",data);
2551 MMG5_SAFE_FREE(data);
2552
2553 if ( ! mesh->adja ) {
2554 if ( !MMG3D_hashTetra(mesh,1) ) {
2555 printf("\n ## Error: %s: unable to compute triangle adjacencies\n.",__func__);
2556 return 0;
2557 }
2558 }
2559
2560 ne = 0;
2561 for (k=1; k<=mesh->ne; k++) {
2562 pt = &mesh->tetra[k];
2563 if ( !MG_EOK(pt) ) continue;
2564 ne++;
2565 }
2566
2567 /* Save elt number, number of neighbors per elt */
2568 fprintf(inm, "%" MMG5_PRId " %d\n\n",ne,mesh->dim+1);
2569
2570 ne = 0;
2571 for ( k=1; k<=mesh->ne; ++k ) {
2572 pt = &mesh->tetra[k];
2573 if ( MG_EOK(pt) ) {
2574 /* Save elt idx */
2575 fprintf(inm, "%" MMG5_PRId " ",++ne);
2576
2577 /* Save neighbors */
2578 for ( i=1; i<=mesh->dim+1; ++i ) {
2579 /* The triangle conventions is that no neighbors <=> -1 */
2580 idx = ( mesh->adja[4*(k-1)+i] > 0 ) ? mesh->adja[4*(k-1)+i]/4 : -1;
2581 fprintf(inm, "%"MMG5_PRId" ",idx);
2582 }
2583 /* Save bdy marker */
2584 fprintf(inm, "\n");
2585 }
2586 }
2587
2588 fclose(inm);
2589
2590 return 1;
2591}
2592
2593static inline
2595 FILE* inm;
2596 MMG5_pTria pt;
2597 MMG5_int k,i;
2598 char *ptr,*data;
2599
2600 if ( !mesh->nt ) {
2601 return 1;
2602 }
2603
2604 if ( (!filename) || !(*filename) ) {
2606 }
2607 if ( (!filename) || !(*filename) ) {
2608 printf("\n ## Error: %s: unable to save a file without a valid filename\n.",
2609 __func__);
2610 return 0;
2611 }
2612
2613 /* Name of file */
2614 MMG5_SAFE_CALLOC(data,strlen(filename)+6,char,return 0);
2615 strcpy(data,filename);
2616 ptr = strstr(data,".node");
2617 if ( ptr ) {
2618 *ptr = '\0';
2619 }
2620
2621 /* Add .node ext */
2622 strcat(data,".face");
2623 if( !(inm = fopen(data,"wb")) ) {
2624 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
2625 MMG5_SAFE_FREE(data);
2626 return 0;
2627 }
2628
2629 fprintf(stdout," %%%% %s OPENED\n",data);
2630 MMG5_SAFE_FREE(data);
2631
2632 /* Save tria number, 1 bdy marker */
2633 fprintf(inm, "%" MMG5_PRId " %d\n\n",mesh->nt,1);
2634
2635 for ( k=1; k<=mesh->nt; ++k ) {
2636 pt = &mesh->tria[k];
2637 if ( MG_EOK(pt) ) {
2638 /* Save elt idx */
2639 fprintf(inm, "%" MMG5_PRId " ",k);
2640
2641 /* Save connectivity */
2642 for ( i=0; i<mesh->dim; ++i ) {
2643 fprintf(inm, "%" MMG5_PRId " ",mesh->point[pt->v[i]].tmp);
2644 }
2645
2646 /* Save bdy marker */
2647 fprintf(inm, "%" MMG5_PRId "\n",pt->ref);
2648 }
2649 }
2650 fprintf(stdout," NUMBER OF TRIANGLES %8" MMG5_PRId "\n",mesh->nt);
2651
2652 fclose(inm);
2653
2654 return 1;
2655}
2656
2657
2659
2660 if ( !MMG5_saveNode(mesh,filename) ) {
2661 return 0;
2662 }
2663
2664 if ( !MMG3D_saveEle(mesh,filename) ) {
2665 return 0;
2666 }
2667
2668 if ( !MMG3D_saveFace(mesh,filename) ) {
2669 return 0;
2670 }
2671
2672 if ( !MMG5_saveEdge(mesh,filename,".edge") ) {
2673 return 0;
2674 }
2675
2676 if ( !MMG3D_saveNeigh(mesh,filename) ) {
2677 return 0;
2678 }
2679
2680 return 1;
2681}
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 MMG3D_Set_inputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solin)
Set the name of input solution file.
int MMG3D_Set_solSize(MMG5_pMesh mesh, MMG5_pSol sol, int typEntity, MMG5_int np, int typSol)
Initialize a solution field.
int ier
tmp[*strlen0]
MMG5_pMesh MMG5_pSol * sol
MMG5_pMesh * mesh
MMG5_pMesh char * filename
int MMG3D_hashTetra(MMG5_pMesh mesh, int pack)
Create array of adjacency.
Definition: hash_3d.c:122
int MMG5_readDoubleSol3D(MMG5_pSol sol, FILE *inm, int bin, int iswp, MMG5_int pos)
Definition: inout.c:2271
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_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
void MMG5_writeDoubleSol3D(MMG5_pMesh mesh, MMG5_pSol sol, FILE *inm, int bin, MMG5_int pos, int metricData)
Definition: inout.c:2320
int MMG5_saveSolAtTetrahedraHeader(MMG5_pMesh mesh, FILE *inm, int ver, int bin, MMG5_int *bpos, int nsols, int nsolsAtTetra, int *entities, int *type, int *size)
Definition: inout.c:2595
int MMG5_readFloatSol3D(MMG5_pSol sol, FILE *inm, int bin, int iswp, int pos)
Definition: inout.c:2222
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
int MMG3D_saveMshMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Save a mesh and data in MSH format, ascii or binary depending on the filename extension.
Definition: inout_3d.c:2138
int MMG3D_openMesh(int imprim, const char *filename, FILE **inm, int *bin, char *modeASCII, char *modeBIN)
Definition: inout_3d.c:53
int MMG3D_saveMesh(MMG5_pMesh mesh, const char *filename)
Save a mesh in .mesh/.meshb format.
Definition: inout_3d.c:1273
int MMG3D_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_3d.c:2062
int MMG3D_loadSol(MMG5_pMesh mesh, MMG5_pSol met, const char *filename)
Load a metric field (or other solution).
Definition: inout_3d.c:2143
int MMG3D_loadMesh(MMG5_pMesh mesh, const char *filename)
Load a mesh (in .mesh/.mesb format) from file.
Definition: inout_3d.c:1049
int MMG3D_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_3d.c:2216
int MMG3D_loadMshMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Load a mesh and possibly a solution in .msh format from file.
Definition: inout_3d.c:1063
int MMG3D_saveTetgenMesh(MMG5_pMesh mesh, const char *filename)
Save data in Tetgen's Triangle format.
Definition: inout_3d.c:2658
static int MMG3D_saveEle(MMG5_pMesh mesh, const char *filename)
Definition: inout_3d.c:2442
static int MMG3D_saveNeigh(MMG5_pMesh mesh, const char *filename)
Definition: inout_3d.c:2514
int MMG3D_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_3d.c:1132
static int MMG3D_saveFace(MMG5_pMesh mesh, const char *filename)
Definition: inout_3d.c:2594
int MMG3D_saveSol(MMG5_pMesh mesh, MMG5_pSol met, const char *filename)
Write isotropic or anisotropic metric.
Definition: inout_3d.c:2314
int MMG3D_saveMshMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Save a mesh in MSH format, ascii or binary depending on the filename extension.
Definition: inout_3d.c:2133
int MMG3D_loadMesh_opened(MMG5_pMesh mesh, FILE *inm, int bin)
Definition: inout_3d.c:113
int MMG3D_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_3d.c:1190
int MMG3D_saveAllSols(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Save 1 or more solutions in medit solution file format.
Definition: inout_3d.c:2352
int MMG3D_loadVtkMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename)
Load a mesh and possibly a solution from a file in VTK format.
int MMG3D_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Save a mesh and optionally one solution in VTK format.
int MMG3D_loadVtuMesh(MMG5_pMesh mesh, MMG5_pSol met, MMG5_pSol sol, const char *filename)
Load a mesh and possibly a solution in VTU (VTK) format from file.
Definition: inoutcpp_3d.cpp:74
int MMG3D_saveVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Save a mesh and optionally one data field in VTU format.
API headers and documentation for the mmg3d library, for volumetric meshes in 3D.
int MMG3D_zaldy(MMG5_pMesh mesh)
Definition: zaldy_3d.c:346
@ 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_Tetrahedron
Definition: libmmgtypes.h:233
#define MG_REQ
#define MG_GEO
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
#define MG_EOK(pt)
#define MG_NUL
#define MG_PARBDY
#define MMG5_SD
#define MMG5_ADD_MEM(mesh, size, message, law)
#define MG_EDG_OR_NOM(tag)
#define MG_SIN(tag)
#define MMG_FREAD(ptr, size, count, stream)
#define MMG5_SW
#define MG_VOK(ppt)
#define MG_CRN
double MMG5_orvol(MMG5_pPoint point, MMG5_int *v)
Definition: tools.c:951
#define MG_BDY
#define MMG5_SAFE_MALLOC(ptr, size, type, law)
#define MG_NOM
#define MMG5_FILESTR_LGTH
#define MMG5_EPSD2
#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
int8_t iso
Definition: libmmgtypes.h:541
MMG mesh structure.
Definition: libmmgtypes.h:613
MMG5_int ntmax
Definition: libmmgtypes.h:620
MMG5_pQuad quadra
Definition: libmmgtypes.h:656
MMG5_int nc1
Definition: libmmgtypes.h:623
MMG5_Info info
Definition: libmmgtypes.h:659
MMG5_int xt
Definition: libmmgtypes.h:628
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_pPrism prism
Definition: libmmgtypes.h:653
MMG5_int nquad
Definition: libmmgtypes.h:621
MMG5_int nemax
Definition: libmmgtypes.h:620
MMG5_int npmax
Definition: libmmgtypes.h:620
MMG5_pxPoint xpoint
Definition: libmmgtypes.h:650
MMG5_int xp
Definition: libmmgtypes.h:628
MMG5_int nei
Definition: libmmgtypes.h:620
MMG5_pTetra tetra
Definition: libmmgtypes.h:651
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 xp
Definition: libmmgtypes.h:285
MMG5_int ref
Definition: libmmgtypes.h:284
MMG5_int flag
Definition: libmmgtypes.h:288
Structure to store prsim of a MMG mesh.
Definition: libmmgtypes.h:469
MMG5_int v[6]
Definition: libmmgtypes.h:470
MMG5_int ref
Definition: libmmgtypes.h:471
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 tetrahedra of an MMG mesh.
Definition: libmmgtypes.h:407
MMG5_int v[4]
Definition: libmmgtypes.h:409
MMG5_int ref
Definition: libmmgtypes.h:410
uint16_t tag
Definition: libmmgtypes.h:417
Structure to store triangles of a MMG mesh.
Definition: libmmgtypes.h:338
MMG5_int ref
Definition: libmmgtypes.h:341
uint16_t tag[3]
Definition: libmmgtypes.h:348
MMG5_int v[3]
Definition: libmmgtypes.h:340
Structure to store surface vertices of an MMG mesh.
Definition: libmmgtypes.h:300
double n1[3]
Definition: libmmgtypes.h:301