Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
inout_s.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
37#include "libmmgs.h"
38#include "libmmgs_private.h"
39#include <math.h>
40
42 FILE *inm;
43 MMG5_pTria pt1,pt2;
44 MMG5_pPoint ppt;
45 double *norm,*n,dd;
46 float fc;
47 long posnp,posnt,posne,posncor,posnq,posned,posnr;
48 long posntreq,posnpreq,posnormal,posnc1;
49 MMG5_int k,ia,nq,nri,ip,idn,ng,npreq,nref;
50 int i;
51 MMG5_int ncor,nedreq,ntreq,posnedreq,bpos;
52 int binch,bin,iswp,bdim;
53 MMG5_int na;
54 char *ptr,*data;
55 char chaine[MMG5_FILESTR_LGTH],strskip[MMG5_FILESTR_LGTH];
56
57 posnp = posnt = posne = posncor = posnq = 0;
58 posned = posnr = posnpreq = posntreq = posnc1 = npreq = 0;
59 posnedreq = posnormal = 0;
60 ncor = nri = ng = nedreq = nq = ntreq = 0;
61 bin = 0;
62 iswp = 0;
63 mesh->np = mesh->nt = mesh->nti = mesh->npi = 0;
64
65 nref = 0;
66
67 MMG5_SAFE_CALLOC(data,strlen(filename)+7,char,return -1);
68
69 strcpy(data,filename);
70 ptr = strstr(data,".mesh");
71 if ( !ptr ) {
72 /* data contains the filename without extension */
73 strcat(data,".meshb");
74 if( !(inm = fopen(data,"rb")) ) {
75 /* our file is not a .meshb file, try with .mesh ext */
76 ptr = strstr(data,".mesh");
77 *ptr = '\0';
78 strcat(data,".mesh");
79 if( !(inm = fopen(data,"rb")) ) {
80 MMG5_SAFE_FREE(data);
81 return 0;
82 }
83 }
84 else bin = 1;
85 }
86 else {
87 ptr = strstr(data,".meshb");
88 if ( ptr ) bin = 1;
89 if( !(inm = fopen(data,"rb")) ) {
90 MMG5_SAFE_FREE(data);
91 return 0;
92 }
93 }
94
95 if ( mesh->info.imprim >= 0 )
96 fprintf(stdout," %%%% %s OPENED\n",data);
97 MMG5_SAFE_FREE(data);
98
99 if (!bin) {
100 strcpy(chaine,"D");
101 while(fscanf(inm,"%127s",&chaine[0])!=EOF && strncmp(chaine,"End",strlen("End")) ) {
102 if ( chaine[0] == '#' ) {
103 fgets(strskip,MMG5_FILESTR_LGTH,inm);
104 continue;
105 }
106 if(!strncmp(chaine,"MeshVersionFormatted",strlen("MeshVersionFormatted"))) {
107 MMG_FSCANF(inm,"%d",&mesh->ver);
108 continue;
109 } else if(!strncmp(chaine,"Dimension",strlen("Dimension"))) {
110 MMG_FSCANF(inm,"%d",&mesh->dim);
111 if(mesh->dim!=3) {
112 fprintf(stderr,"BAD DIMENSION : %d\n",mesh->dim);
113 return -1;
114 }
115 continue;
116 } else if(!strncmp(chaine,"Vertices",strlen("Vertices"))) {
117 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->npi);
118 posnp = ftell(inm);
119 continue;
120 } else if(!strncmp(chaine,"RequiredVertices",strlen("RequiredVertices"))) {
121 MMG_FSCANF(inm,"%" MMG5_PRId "",&npreq);
122 posnpreq = ftell(inm);
123 continue;
124 } else if(!strncmp(chaine,"Triangles",strlen("Triangles"))) {
125 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->nti);
126 posnt = ftell(inm);
127 continue;
128 }
129 else if(!strncmp(chaine,"RequiredTriangles",strlen("RequiredTriangles"))) {
130 MMG_FSCANF(inm,"%" MMG5_PRId "",&ntreq);
131 posntreq = ftell(inm);
132 continue;
133 } else if(!strncmp(chaine,"Quadrilaterals",strlen("Quadrilaterals"))) {
134 MMG_FSCANF(inm,"%" MMG5_PRId "",&nq);
135 posnq = ftell(inm);
136 continue;
137 } else if(!strncmp(chaine,"Corners",strlen("Corners"))) {
138 MMG_FSCANF(inm,"%" MMG5_PRId "",&ncor);
139 posncor = ftell(inm);
140 continue;
141 } else if(!strncmp(chaine,"Edges",strlen("Edges"))) {
142 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->na);
143 posned = ftell(inm);
144 continue;
145 } else if(!strncmp(chaine,"RequiredEdges",strlen("RequiredEdges"))) {
146 MMG_FSCANF(inm,"%" MMG5_PRId "",&nedreq);
147 posnedreq = ftell(inm);
148 continue;
149 } else if(!strncmp(chaine,"Ridges",strlen("Ridges"))) {
150 MMG_FSCANF(inm,"%" MMG5_PRId "",&nri);
151 posnr = ftell(inm);
152 continue;
153 } else if(!ng && !strncmp(chaine,"Normals",strlen("Normals"))) {
154 MMG_FSCANF(inm,"%" MMG5_PRId "",&ng);
155 posnormal = ftell(inm);
156 continue;
157 } else if(!strncmp(chaine,"NormalAtVertices",strlen("NormalAtVertices"))) {
158 MMG_FSCANF(inm,"%" MMG5_PRId "",&mesh->nc1);
159 posnc1 = ftell(inm);
160 continue;
161 }
162 }
163 } else { //binary file
164 bdim = 0;
165 MMG_FREAD(&mesh->ver,MMG5_SW,1,inm);
166 iswp=0;
167 if(mesh->ver==16777216)
168 iswp=1;
169 else if(mesh->ver!=1) {
170 fprintf(stdout,"BAD FILE ENCODING\n");
171 }
172 MMG_FREAD(&mesh->ver,MMG5_SW,1,inm);
173 if(iswp) mesh->ver = MMG5_swapbin(mesh->ver);
174 while(fread(&binch,MMG5_SW,1,inm)!=0 && binch!=54 ) {
175 if(iswp) binch=MMG5_swapbin(binch);
176 if(binch==54) break;
177 if(!bdim && binch==3) { //Dimension
178 MMG_FREAD(&bdim,MMG5_SW,1,inm); //NulPos=>20
179 if(iswp) bdim=MMG5_swapbin(bdim);
180 MMG_FREAD(&bdim,MMG5_SW,1,inm);
181 if(iswp) bdim=MMG5_swapbin(bdim);
182 mesh->dim = bdim;
183 if(bdim!=3) {
184 fprintf(stderr,"BAD MESH DIMENSION : %d\n",mesh->dim);
185 return -1;
186 }
187 continue;
188 } else if(!mesh->npi && binch==4) { //Vertices
189 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
190 if(iswp) bpos=MMG5_swapbin(bpos);
191 MMG_FREAD(&mesh->npi,MMG5_SW,1,inm);
192 if(iswp) mesh->npi=MMG5_SWAPBIN(mesh->npi);
193 posnp = ftell(inm);
194 rewind(inm);
195 fseek(inm,bpos,SEEK_SET);
196 continue;
197 } else if(binch==15) { //RequiredVertices
198 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
199 if(iswp) bpos=MMG5_swapbin(bpos);
200 MMG_FREAD(&npreq,MMG5_SW,1,inm);
201 if(iswp) npreq=MMG5_SWAPBIN(npreq);
202 posnpreq = ftell(inm);
203 rewind(inm);
204 fseek(inm,bpos,SEEK_SET);
205 continue;
206 } else if(!mesh->nti && binch==6) {//Triangles
207 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
208 if(iswp) bpos=MMG5_swapbin(bpos);
209 MMG_FREAD(&mesh->nti,MMG5_SW,1,inm);
210 if(iswp) mesh->nti=MMG5_SWAPBIN(mesh->nti);
211 posnt = ftell(inm);
212 rewind(inm);
213 fseek(inm,bpos,SEEK_SET);
214 continue;
215 } else if(binch==17) { //RequiredTriangles
216 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
217 if(iswp) bpos=MMG5_swapbin(bpos);
218 MMG_FREAD(&ntreq,MMG5_SW,1,inm);
219 if(iswp) ntreq=MMG5_SWAPBIN(ntreq);
220 posntreq = ftell(inm);
221 rewind(inm);
222 fseek(inm,bpos,SEEK_SET);
223 continue;
224 } else if(binch==7) {//Quadrilaterals
225 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
226 if(iswp) bpos=MMG5_swapbin(bpos);
227 MMG_FREAD(&nq,MMG5_SW,1,inm);
228 if(iswp) nq=MMG5_SWAPBIN(nq);
229 posnq = ftell(inm);
230 rewind(inm);
231 fseek(inm,bpos,SEEK_SET);
232 continue;
233 } else if(!ncor && binch==13) { //Corners
234 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
235 if(iswp) bpos=MMG5_swapbin(bpos);
236 MMG_FREAD(&ncor,MMG5_SW,1,inm);
237 if(iswp) ncor=MMG5_SWAPBIN(ncor);
238 posncor = ftell(inm);
239 rewind(inm);
240 fseek(inm,bpos,SEEK_SET);
241 continue;
242 } else if(!mesh->na && binch==5) { //Edges
243 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
244 if(iswp) bpos=MMG5_swapbin(bpos);
245 MMG_FREAD(&mesh->na,MMG5_SW,1,inm);
246 if(iswp) mesh->na=MMG5_SWAPBIN(mesh->na);
247 posned = ftell(inm);
248 rewind(inm);
249 fseek(inm,bpos,SEEK_SET);
250 continue;
251 } else if(binch==16) { //RequiredEdges
252 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
253 if(iswp) bpos=MMG5_swapbin(bpos);
254 MMG_FREAD(&nedreq,MMG5_SW,1,inm);
255 if(iswp) nedreq=MMG5_SWAPBIN(nedreq);
256 posnedreq = ftell(inm);
257 rewind(inm);
258 fseek(inm,bpos,SEEK_SET);
259 continue;
260 } else if(binch==14) { //Ridges
261 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
262 if(iswp) bpos=MMG5_swapbin(bpos);
263 MMG_FREAD(&nri,MMG5_SW,1,inm);
264 if(iswp) nri=MMG5_SWAPBIN(nri);
265 posnr = ftell(inm);
266 rewind(inm);
267 fseek(inm,bpos,SEEK_SET);
268 continue;
269 } else if(!ng && binch==60) { //Normals
270 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
271 if(iswp) bpos=MMG5_swapbin(bpos);
272 MMG_FREAD(&ng,MMG5_SW,1,inm);
273 if(iswp) ng=MMG5_SWAPBIN(ng);
274 posnormal = ftell(inm);
275 rewind(inm);
276 fseek(inm,bpos,SEEK_SET);
277 continue;
278 } else if(binch==20) { //NormalAtVertices
279 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
280 if(iswp) bpos=MMG5_swapbin(bpos);
281 MMG_FREAD(&mesh->nc1,MMG5_SW,1,inm);
282 if(iswp) mesh->nc1=MMG5_SWAPBIN(mesh->nc1);
283 posnc1 = ftell(inm);
284 rewind(inm);
285 fseek(inm,bpos,SEEK_SET);
286 continue;
287 } else {
288 MMG_FREAD(&bpos,MMG5_SW,1,inm); //NulPos
289 if(iswp) bpos=MMG5_swapbin(bpos);
290 rewind(inm);
291 fseek(inm,bpos,SEEK_SET);
292 }
293 }
294 }
295
296 if ( !mesh->npi || !mesh->nti ) {
297 fprintf(stdout," ** MISSING DATA\n");
298 return -1;
299 }
300 mesh->np = mesh->npi;
301 mesh->nt = mesh->nti + 2*nq;
302
303 /* mem alloc */
304 if ( !MMGS_zaldy(mesh) ) return -1;
305
306 /* read vertices */
307
308 rewind(inm);
309 fseek(inm,posnp,SEEK_SET);
310 for (k=1; k<=mesh->np; k++) {
311 ppt = &mesh->point[k];
312 if (mesh->ver < 2) { /*float*/
313 if (!bin) {
314 for (i=0 ; i<3 ; i++) {
315 MMG_FSCANF(inm,"%f",&fc);
316 ppt->c[i] = (double) fc;
317 }
318 MMG_FSCANF(inm,"%" MMG5_PRId "",&ppt->ref);
319 } else {
320 for (i=0 ; i<3 ; i++) {
321 MMG_FREAD(&fc,MMG5_SW,1,inm);
322 if(iswp) fc=MMG5_swapf(fc);
323 ppt->c[i] = (double) fc;
324 }
325 MMG_FREAD(&ppt->ref,MMG5_SW,1,inm);
326 if(iswp) ppt->ref=MMG5_SWAPBIN(ppt->ref);
327 }
328 } else {
329 if (!bin) {
330 MMG_FSCANF(inm,"%lf %lf %lf %" MMG5_PRId "",&ppt->c[0],&ppt->c[1],&ppt->c[2],&ppt->ref);
331 }
332 else {
333 for (i=0 ; i<3 ; i++) {
334 MMG_FREAD(&ppt->c[i],MMG5_SD,1,inm);
335 if(iswp) ppt->c[i]=MMG5_swapd(ppt->c[i]);
336 }
337 MMG_FREAD(&ppt->ref,MMG5_SW,1,inm);
338 if(iswp) ppt->ref=MMG5_SWAPBIN(ppt->ref);
339 }
340 }
341 if ( ppt->ref < 0 ) {
342 ppt->ref = -ppt->ref;
343 ++nref;
344 }
345 ppt->tag = MG_NUL;
346 }
347
348 /* read triangles and set seed */
349 if ( mesh->nti ) {
350 rewind(inm);
351 fseek(inm,posnt,SEEK_SET);
352 for (k=1; k<=mesh->nti; k++) {
353 pt1 = &mesh->tria[k];
354 if (!bin) {
355 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&pt1->v[0],&pt1->v[1],&pt1->v[2],&pt1->ref);
356 }
357 else {
358 for (i=0 ; i<3 ; i++) {
359 MMG_FREAD(&pt1->v[i],MMG5_SW,1,inm);
360 if(iswp) pt1->v[i]=MMG5_SWAPBIN(pt1->v[i]);
361 }
362 MMG_FREAD(&pt1->ref,MMG5_SW,1,inm);
363 if(iswp) pt1->ref=MMG5_SWAPBIN(pt1->ref);
364 }
365 for (i=0; i<3; i++) {
366 ppt = &mesh->point[pt1->v[i]];
367 ppt->tag &= ~MG_NUL;
368 }
369 }
370 /* get required triangles */
371 if(ntreq) {
372 rewind(inm);
373 fseek(inm,posntreq,SEEK_SET);
374 for (k=1; k<=ntreq; k++) {
375 if(!bin) {
376 MMG_FSCANF(inm,"%d",&i);
377 }
378 else {
379 MMG_FREAD(&i,MMG5_SW,1,inm);
380 if(iswp) i=MMG5_swapbin(i);
381 }
382 if ( i>mesh->nti ) {
383 fprintf(stderr,"\n ## Warning: %s: required triangle number %8d"
384 " ignored.\n",__func__,i);
385 } else {
386 pt1 = &mesh->tria[i];
387 pt1->tag[0] |= MG_REQ;
388 pt1->tag[1] |= MG_REQ;
389 pt1->tag[2] |= MG_REQ;
390 }
391 }
392 }
393 }
394
395 /* read quads: automatic conversion into tria */
396 if ( nq > 0 ) {
397 rewind(inm);
398 fseek(inm,posnq,SEEK_SET);
399
400 printf(" ## Warning: %s: quadrangles automatically converted into"
401 " triangles\n.",__func__);
402
403 for (k=1; k<=nq; k++) {
404 mesh->nti++;
405 pt1 = &mesh->tria[mesh->nti];
406 mesh->nti++;
407 pt2 = &mesh->tria[mesh->nti];
408
409 if (!bin) {
410 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&pt1->v[0],&pt1->v[1],&pt1->v[2],&pt2->v[2],&pt1->ref);
411 }
412 else {
413 for (i=0 ; i<3 ; i++) {
414 MMG_FREAD(&pt1->v[i],MMG5_SW,1,inm);
415 if(iswp) pt1->v[i]=MMG5_SWAPBIN(pt1->v[i]);
416 }
417 MMG_FREAD(&pt2->v[2],MMG5_SW,1,inm);
418 if(iswp) pt2->v[2]=MMG5_SWAPBIN(pt2->v[2]);
419 MMG_FREAD(&pt1->ref,MMG5_SW,1,inm);
420 if(iswp) pt1->ref=MMG5_SWAPBIN(pt1->ref);
421 }
422 if ( pt1->ref < 0 ) {
423 pt1->ref = -pt1->ref;
424 nref += 2;
425 }
426
427 pt2->v[0] = pt1->v[0];
428 pt2->v[1] = pt1->v[2];
429 pt2->ref = pt1->ref;
430 for (i=0; i<3; i++) {
431 ppt = &mesh->point[pt1->v[i]];
432 ppt->tag &= ~MG_NUL;
433 ppt = &mesh->point[pt2->v[i]];
434 ppt->tag &= ~MG_NUL;
435 }
436 }
437 mesh->nt = mesh->nti;
438 }
439
440 if(ncor) {
441 rewind(inm);
442 fseek(inm,posncor,SEEK_SET);
443 for (k=1; k<=ncor; k++) {
444 if(!bin) {
445 MMG_FSCANF(inm,"%d",&i);
446 }
447 else {
448 MMG_FREAD(&i,MMG5_SW,1,inm);
449 if(iswp) i=MMG5_swapbin(i);
450 }
451 if(i>mesh->np) {
452 fprintf(stderr,"\n ## Warning: %s: corner number %8d ignored.\n",
453 __func__,i);
454 } else {
455 ppt = &mesh->point[i];
456 ppt->tag |= MG_CRN;
457 }
458 }
459 }
460
461 /* read required vertices */
462 if(npreq) {
463 rewind(inm);
464 fseek(inm,posnpreq,SEEK_SET);
465 for (k=1; k<=npreq; k++) {
466 if(!bin) {
467 MMG_FSCANF(inm,"%d",&i);
468 }
469 else {
470 MMG_FREAD(&i,MMG5_SW,1,inm);
471 if(iswp) i=MMG5_swapbin(i);
472 }
473 if(i>mesh->np) {
474 fprintf(stderr,"\n ## Warning: %s: required Vertices number %8d ignored\n",
475 __func__,i);
476 } else {
477 ppt = &mesh->point[i];
478 ppt->tag |= MG_REQ;
479 ppt->tag &= ~MG_NUL;
480 }
481 }
482 }
483
484 /* Read mesh edges */
485 na = mesh->na;
486 if ( mesh->na ) {
487 rewind(inm);
488 fseek(inm,posned,SEEK_SET);
489
490 for (k=1; k<=mesh->na; k++) {
491 if (!bin) {
492 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&mesh->edge[k].a,&mesh->edge[k].b,&mesh->edge[k].ref);
493 }
494 else {
495 MMG_FREAD(&mesh->edge[k].a,MMG5_SW,1,inm);
496 if(iswp) mesh->edge[k].a=MMG5_SWAPBIN(mesh->edge[k].a);
497 MMG_FREAD(&mesh->edge[k].b,MMG5_SW,1,inm);
498 if(iswp) mesh->edge[k].b=MMG5_SWAPBIN(mesh->edge[k].b);
499 MMG_FREAD(&mesh->edge[k].ref,MMG5_SW,1,inm);
500 if(iswp) mesh->edge[k].ref=MMG5_SWAPBIN(mesh->edge[k].ref);
501 }
502 if ( mesh->edge[k].ref < 0 ) {
503 mesh->edge[k].ref = -mesh->edge[k].ref;
504 ++nref;
505 }
506 mesh->edge[k].tag |= MG_REF;
507 mesh->point[mesh->edge[k].a].tag |= MG_REF;
508 mesh->point[mesh->edge[k].b].tag |= MG_REF;
509 }
510
511 /* get ridges */
512 if ( nri ) {
513 rewind(inm);
514 fseek(inm,posnr,SEEK_SET);
515 for (k=1; k<=nri; k++) {
516 if(!bin) {
517 MMG_FSCANF(inm,"%" MMG5_PRId "",&ia);
518 }
519 else {
520 MMG_FREAD(&ia,MMG5_SW,1,inm);
521 if(iswp) ia=MMG5_SWAPBIN(ia);
522 }
523 if ( (ia>na) || (ia<0) ) {
524 fprintf(stderr,"\n ## Warning: %s: ridge number %8" MMG5_PRId " ignored.\n",
525 __func__,ia);
526 }
527 else {
528 mesh->edge[ia].tag |= MG_GEO;
529 }
530 }
531 }
532
533 if ( nedreq ) {
534 rewind(inm);
535 fseek(inm,posnedreq,SEEK_SET);
536 for (k=1; k<=nedreq; k++) {
537 if(!bin) {
538 MMG_FSCANF(inm,"%" MMG5_PRId "",&ia);
539 }
540 else {
541 MMG_FREAD(&ia,MMG5_SW,1,inm);
542 if(iswp) ia=MMG5_SWAPBIN(ia);
543 }
544 if ( (ia>na) || (ia<0) ) {
545 fprintf(stderr,"\n ## Warning: %s: required edge number %8" MMG5_PRId " ignored\n",
546 __func__,ia);
547 }
548 else {
549 mesh->edge[ia].tag |= MG_REQ;
550 }
551 }
552 }
553 }
554
555 /* read geometric entities */
556 if ( mesh->nc1 && !ng ) {
557 fprintf(stderr,"\n ## Warning: %s: your mesh don't contains Normals but contains"
558 " NormalAtVertices. The NormalAtVertices are deleted. \n",__func__);
559 mesh->nc1 = 0;
560 }
561
562 if ( ng > 0 ) {
563 MMG5_SAFE_CALLOC(norm,3*ng+1,double,return -1);
564
565 rewind(inm);
566 fseek(inm,posnormal,SEEK_SET);
567 for (k=1; k<=ng; k++) {
568 n = &norm[3*(k-1)+1];
569 if ( mesh->ver == 1 ) {
570 if (!bin) {
571 for (i=0 ; i<3 ; i++) {
572 MMG_FSCANF(inm,"%f",&fc);
573 n[i] = (double) fc;
574 }
575 } else {
576 for (i=0 ; i<3 ; i++) {
577 MMG_FREAD(&fc,MMG5_SW,1,inm);
578 if(iswp) fc=MMG5_swapf(fc);
579 n[i] = (double) fc;
580 }
581 }
582 }
583 else {
584 if (!bin) {
585 MMG_FSCANF(inm,"%lf %lf %lf",&n[0],&n[1],&n[2]);
586 }
587 else {
588 for (i=0 ; i<3 ; i++) {
589 MMG_FREAD(&n[i],MMG5_SD,1,inm);
590 if(iswp) n[i]=MMG5_swapd(n[i]);
591 }
592 }
593 }
594 dd = n[0]*n[0] + n[1]*n[1] + n[2]*n[2];
595 if ( dd > MMG5_EPSD2 ) {
596 dd = 1.0 / sqrt(dd);
597 n[0] *= dd;
598 n[1] *= dd;
599 n[2] *= dd;
600 }
601 }
602
603 rewind(inm);
604 fseek(inm,posnc1,SEEK_SET);
605
606 for (k=1; k<=mesh->nc1; k++) {
607 if (!bin) {
608 MMG_FSCANF(inm,"%" MMG5_PRId " %" MMG5_PRId "",&ip,&idn);
609 }
610 else {
611 MMG_FREAD(&ip,MMG5_SW,1,inm);
612 if(iswp) ip=MMG5_SWAPBIN(ip);
613 MMG_FREAD(&idn,MMG5_SW,1,inm);
614 if(iswp) idn=MMG5_SWAPBIN(idn);
615 }
616 if ( idn > 0 && ip < mesh->np+1 )
617 memcpy(&mesh->point[ip].n,&norm[3*(idn-1)+1],3*sizeof(double));
618 }
619 MMG5_SAFE_FREE(norm);
620 }
621
622 if ( nref ) {
623 fprintf(stdout,"\n $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n");
624 fprintf(stdout," WARNING : %" MMG5_PRId " entities with unexpected refs (ref< 0).\n",nref);
625 fprintf(stdout," We take their absolute values.\n");
626 fprintf(stdout," $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n\n");
627 }
628
629 if ( abs(mesh->info.imprim) > 4 ) {
630 fprintf(stdout," NUMBER OF VERTICES %8" MMG5_PRId " / %8" MMG5_PRId " CORNERS/REQ. %" MMG5_PRId " / %" MMG5_PRId "\n",mesh->npi,mesh->npmax,ncor,npreq);
631 fprintf(stdout," NUMBER OF TRIANGLES %8" MMG5_PRId " / %8" MMG5_PRId "\n",mesh->nti,mesh->ntmax);
632
633 if ( mesh->na )
634 fprintf(stdout," NUMBER OF EDGES %8" MMG5_PRId " RIDGES %6" MMG5_PRId "\n",mesh->na,nri);
635 }
636 fclose(inm);
637 return 1;
638}
639
641 FILE* inm;
642 long posNodes,posElts,*posNodeData;
643 int ier,bin,iswp,nsols;
644 MMG5_int nelts;
645
646 mesh->dim = 3;
647
649 &posNodes,&posElts,&posNodeData,
650 &bin,&iswp,&nelts,&nsols);
651 if ( ier < 1 ) return (ier);
652
653 if ( nsols > 1 ) {
654 fprintf(stderr,"Error: SEVERAL SOLUTIONS FOUND (%d)\n",nsols);
655 fclose(inm);
656 MMG5_SAFE_FREE(posNodeData);
657 return -1;
658 }
659
660 if ( !MMGS_zaldy(mesh) ) {
661 fclose(inm);
662 MMG5_SAFE_FREE(posNodeData);
663 return -1;
664 }
665
666 mesh->ne = mesh->nprism = 0;
667
668 if ( !mesh->nt ) {
669 fprintf(stderr," ** MISSING DATA.\n");
670 fprintf(stderr," Check that your mesh contains triangles.\n");
671 fprintf(stderr," Exit program.\n");
672 fclose(inm);
673 MMG5_SAFE_FREE(posNodeData);
674 return -1;
675 }
676
677
678 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt ) {
679 fclose(inm);
680 MMG5_SAFE_FREE(posNodeData);
681 return -1;
682 }
683
685 posNodes,posElts,posNodeData,
686 bin,iswp,nelts,nsols);
687 MMG5_SAFE_FREE(posNodeData);
688 if ( ier < 1 ) {
689 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
690 return ier;
691 }
692
693 if ( sol ) {
694 /* Check the metric type */
695 ier = MMG5_chkMetricType(mesh,&sol->type,&sol->entities,inm);
696 }
697
698 return ier;
699}
700
702 FILE* inm;
703 long posNodes,posElts,*posNodeData;
704 int ier,bin,iswp,nsols;
705 MMG5_int nelts;
706
707 mesh->dim = 3;
708
710 &posNodes,&posElts,&posNodeData,
711 &bin,&iswp,&nelts,&nsols);
712 if ( ier < 1 ) return (ier);
713
714 mesh->nsols = nsols;
715 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
716 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
717 printf(" Exit program.\n"); fclose(inm);
718 MMG5_SAFE_FREE(posNodeData);
719 return -1);
720 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
721
722 if ( !MMGS_zaldy(mesh) ) {
723 fclose(inm);
724 MMG5_SAFE_FREE(posNodeData);
725 return -1;
726 }
727
728 mesh->ne = mesh->nprism = 0;
729
730 if ( !mesh->nt ) {
731 fprintf(stderr," ** MISSING DATA.\n");
732 fprintf(stderr," Check that your mesh contains triangles.\n");
733 fprintf(stderr," Exit program.\n");
734 fclose(inm);
735 MMG5_SAFE_FREE(posNodeData);
736 return -1;
737 }
738
739
740 if (mesh->npmax < mesh->np || mesh->ntmax < mesh->nt ) {
741 fclose(inm);
742 MMG5_SAFE_FREE(posNodeData);
743 return -1;
744 }
745
747 posNodes,posElts,posNodeData,
748 bin,iswp,nelts,nsols);
749
750 if ( ier < 1 ) {
751 fprintf(stderr," ** ERROR WHEN PARSING THE INPUT FILE\n");
752 }
753
754 MMG5_SAFE_FREE(posNodeData);
755
756 return ier;
757}
758
760 int ier=0;
761 const char *filenameptr,*solnameptr;
762 char *tmp,*soltmp;
763
764 if ( filename && strlen(filename) ) {
765 filenameptr = filename;
766 solnameptr = filename;
767 }
768 else if (mesh->namein && strlen(mesh->namein) ) {
769 filenameptr = mesh->namein;
770 if ( sol && strlen(sol->namein) ) {
771 solnameptr = sol->namein;
772 }
773 else {
774 solnameptr = mesh->namein;
775 }
776 }
777 else {
778 fprintf(stderr," ## Error: %s: please provide input file name"
779 " (either in the mesh structure or as function argument).\n",
780 __func__);
781 return -1;
782 }
783
784 MMG5_SAFE_MALLOC(tmp,strlen(filenameptr)+1,char,return -1);
785 strcpy(tmp,filenameptr);
786
787 /* read mesh/sol files */
788 char *ptr = MMG5_Get_filenameExt(tmp);
789 int fmtin = MMG5_Get_format(ptr,MMG5_FMT_MeditASCII);
790
791 switch ( fmtin ) {
792
793 case ( MMG5_FMT_GmshASCII ): case ( MMG5_FMT_GmshBinary ):
795 break;
796
797 case ( MMG5_FMT_VtkVtp ):
799 break;
800
801 case ( MMG5_FMT_VtkVtu ):
803 break;
804
805 case ( MMG5_FMT_VtkVtk ):
807 break;
808
809 case ( MMG5_FMT_MeditASCII ): case ( MMG5_FMT_MeditBinary ):
811 if ( ier < 1 ) { break; }
812
813 /* Facultative metric */
814 if ( sol ) {
815 MMG5_SAFE_MALLOC(soltmp,strlen(solnameptr)+1,char,return -1);
816 strcpy(soltmp,solnameptr);
817
818 if ( MMGS_loadSol(mesh,sol,tmp) == -1) {
819 fprintf(stderr,"\n ## ERROR: WRONG DATA TYPE OR WRONG SOLUTION NUMBER.\n");
820 ier = 0;
821 }
822 MMG5_SAFE_FREE(soltmp);
823 }
824 break;
825
826 default:
827 fprintf(stderr," ** I/O AT FORMAT %s NOT IMPLEMENTED.\n",MMG5_Get_formatName(fmtin) );
828 ier= 0;
829 }
830
832
833 return ier;
834}
835
836
838 FILE *inm;
839 MMG5_pPoint ppt;
840 MMG5_pTria pt;
841 MMG5_pxPoint go;
842 MMG5_int k,np,nt,nc,ng,nn,nr,nre,ntreq,bpos;
843 int bin,binch;
844 // int outm;
845 char *data,*ptr,chaine[MMG5_FILESTR_LGTH];
846
847 mesh->ver = 2;
848
849 bin = 0;
850
851 MMG5_SAFE_CALLOC(data,strlen(filename)+7,char,return 0);
852
853 strcpy(data,filename);
854 ptr = strstr(data,".mesh");
855 if ( !ptr ) {
856 strcat(data,".meshb");
857 if( !(inm = fopen(data,"wb")) ) {
858 ptr = strstr(data,".mesh");
859 *ptr = '\0';
860 strcat(data,".mesh");
861 if( !(inm = fopen(data,"w")) ) {
862 MMG5_SAFE_FREE(data);
863 return 0;
864 }
865 } else {
866 bin = 1;
867 }
868 }
869 else {
870 ptr = strstr(data,".meshb");
871 if( ptr ) {
872 bin = 1;
873 if( !(inm = fopen(data,"wb")) ) {
874 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
875 MMG5_SAFE_FREE(data);
876 return 0;
877 }
878 } else {
879 if( !(inm = fopen(data,"w")) ) {
880 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
881 MMG5_SAFE_FREE(data);
882 return 0;
883 }
884 }
885 }
886 if ( mesh->info.imprim >= 0 )
887 fprintf(stdout," %%%% %s OPENED\n",data);
888 MMG5_SAFE_FREE(data);
889
890 /*entete fichier*/
891 if(!bin) {
892 strcpy(&chaine[0],"MeshVersionFormatted 2\n");
893 fprintf(inm,"%s",chaine);
894 strcpy(&chaine[0],"\n\nDimension 3\n");
895 fprintf(inm,"%s ",chaine);
896 } else {
897 binch = 1; //MeshVersionFormatted
898 fwrite(&binch,MMG5_SW,1,inm);
899 binch = 2; //version
900 fwrite(&binch,MMG5_SW,1,inm);
901 binch = 3; //Dimension
902 fwrite(&binch,MMG5_SW,1,inm);
903 bpos = 20; //Pos
904 fwrite(&bpos,MMG5_SW,1,inm);
905 binch = 3;
906 fwrite(&binch,MMG5_SW,1,inm);
907
908 }
909 /* vertices */
910 np = nc = ng = nn = nre = 0;
911 for (k=1; k<=mesh->np; k++) {
912 ppt = &mesh->point[k];
913 ppt->tmp = 0;
914 if ( MG_VOK(ppt) ) {
915 np++;
916 ppt->tmp = np;
917 if ( ppt->tag & MG_CRN ) nc++;
918 if ( ppt->tag & MG_REQ ) nre++;
919 if ( MG_EDG(ppt->tag) ) ng++;
920 }
921 }
922
923 if(!bin) {
924 strcpy(&chaine[0],"\n\nVertices\n");
925 fprintf(inm,"%s",chaine);
926 fprintf(inm,"%" MMG5_PRId "\n",np);
927 } else {
928 binch = 4; //Vertices
929 fwrite(&binch,MMG5_SW,1,inm);
930 bpos += (3+(1+3*mesh->ver)*np)*MMG5_SW; //NullPos
931 fwrite(&bpos,MMG5_SW,1,inm);
932 fwrite(&np,MMG5_SW,1,inm);
933 }
934 for (k=1; k<=mesh->np; k++) {
935 ppt = &mesh->point[k];
936 if ( MG_VOK(ppt) ) {
937 if(!bin) {
938 fprintf(inm,"%.15lg %.15lg %.15lg %" MMG5_PRId "\n",ppt->c[0],ppt->c[1],ppt->c[2],MMG5_abs(ppt->ref));
939 } else {
940 fwrite((unsigned char*)&ppt->c[0],MMG5_SD,1,inm);
941 fwrite((unsigned char*)&ppt->c[1],MMG5_SD,1,inm);
942 fwrite((unsigned char*)&ppt->c[2],MMG5_SD,1,inm);
943 ppt->ref = MMG5_abs(ppt->ref);
944 fwrite((unsigned char*)&ppt->ref,MMG5_SW,1,inm);
945 }
946 if ( !((ppt->tag & MG_GEO) || ppt->tag & MG_CRN) ) nn++;
947 }
948 }
949
950 nt = 0;
951 for (k=1; k<=mesh->nt; k++) {
952 pt = &mesh->tria[k];
953 if ( !MG_EOK(pt) ) continue;
954 nt++;
955 }
956
957 /* write corners */
958 if ( nc ) {
959 if(!bin) {
960 strcpy(&chaine[0],"\n\nCorners\n");
961 fprintf(inm,"%s",chaine);
962 fprintf(inm,"%" MMG5_PRId "\n",nc);
963 } else {
964 binch = 13; //
965 fwrite(&binch,MMG5_SW,1,inm);
966 bpos += (3+nc)*MMG5_SW; //NullPos
967 fwrite(&bpos,MMG5_SW,1,inm);
968 fwrite(&nc,MMG5_SW,1,inm);
969 }
970 for (k=1; k<=mesh->np; k++) {
971 ppt = &mesh->point[k];
972 if ( MG_VOK(ppt) && ppt->tag & MG_CRN ) {
973 if(!bin) {
974 fprintf(inm,"%" MMG5_PRId "\n",ppt->tmp);
975 } else {
976 fwrite(&ppt->tmp,MMG5_SW,1,inm);
977 }
978 }
979 }
980 }
981 if ( nre ) {
982 if(!bin) {
983 strcpy(&chaine[0],"\n\nRequiredVertices\n");
984 fprintf(inm,"%s",chaine);
985 fprintf(inm,"%" MMG5_PRId "\n",nre);
986 } else {
987 binch = 15; //
988 fwrite(&binch,MMG5_SW,1,inm);
989 bpos += (3+nre)*MMG5_SW; //NullPos
990 fwrite(&bpos,MMG5_SW,1,inm);
991 fwrite(&nre,MMG5_SW,1,inm);
992 }
993 for (k=1; k<=mesh->np; k++) {
994 ppt = &mesh->point[k];
995 if ( MG_VOK(ppt) && ppt->tag & MG_REQ ) {
996 if(!bin) {
997 fprintf(inm,"%" MMG5_PRId "\n",ppt->tmp);
998 } else {
999 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1000 }
1001 }
1002 }
1003 }
1004
1005 /* write edges, ridges */
1006 if ( mesh->na ) {
1007 if(!bin) {
1008 strcpy(&chaine[0],"\n\nEdges\n");
1009 fprintf(inm,"%s",chaine);
1010 fprintf(inm,"%" MMG5_PRId "\n",mesh->na);
1011 } else {
1012 binch = 5; //Edges
1013 fwrite(&binch,MMG5_SW,1,inm);
1014 bpos += (3 + 3*mesh->na)*MMG5_SW;//Pos
1015 fwrite(&bpos,MMG5_SW,1,inm);
1016 fwrite(&mesh->na,MMG5_SW,1,inm);
1017 }
1018 nre = nr = 0;
1019 for (k=1; k<=mesh->na; k++) {
1020 if(!bin) {
1021 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",
1022 mesh->edge[k].a,mesh->edge[k].b,mesh->edge[k].ref);
1023 } else {
1024 fwrite(&mesh->edge[k].a,MMG5_SW,1,inm);
1025 fwrite(&mesh->edge[k].b,MMG5_SW,1,inm);
1026 fwrite(&mesh->edge[k].ref,MMG5_SW,1,inm);
1027 }
1028 if ( mesh->edge[k].tag & MG_REQ ) nre++;
1029 if ( mesh->edge[k].tag & MG_GEO ) nr++;
1030 }
1031 /* ridges */
1032 if ( nr ) {
1033 if(!bin) {
1034 strcpy(&chaine[0],"\n\nRidges\n");
1035 fprintf(inm,"%s",chaine);
1036 fprintf(inm,"%" MMG5_PRId "\n",nr);
1037 } else {
1038 binch = 14; //Ridges
1039 fwrite(&binch,MMG5_SW,1,inm);
1040 bpos += (3 + nr)*MMG5_SW;//Pos
1041 fwrite(&bpos,MMG5_SW,1,inm);
1042 fwrite(&nr,MMG5_SW,1,inm);
1043 }
1044 for (k=1; k<=mesh->na; k++) {
1045 if ( mesh->edge[k].tag & MG_GEO ) {
1046 if(!bin) {
1047 fprintf(inm,"%" MMG5_PRId "\n",k);
1048 } else {
1049 fwrite(&k,MMG5_SW,1,inm);
1050 }
1051 }
1052 }
1053 }
1054 if ( nre ) {
1055 if(!bin) {
1056 strcpy(&chaine[0],"\n\nRequiredEdges\n");
1057 fprintf(inm,"%s",chaine);
1058 fprintf(inm,"%" MMG5_PRId "\n",nre);
1059 } else {
1060 binch = 16; //RequiredEdges
1061 fwrite(&binch,MMG5_SW,1,inm);
1062 bpos += (3 + nre)*MMG5_SW;//Pos
1063 fwrite(&bpos,MMG5_SW,1,inm);
1064 fwrite(&nre,MMG5_SW,1,inm);
1065 }
1066 for (k=1; k<=mesh->na; k++)
1067 if ( mesh->edge[k].tag & MG_REQ ) {
1068 if(!bin) {
1069 fprintf(inm,"%" MMG5_PRId "\n",k);
1070 } else {
1071 fwrite(&k,MMG5_SW,1,inm);
1072 }
1073 }
1074 }
1075 }
1076
1077 /* write triangles */
1078 if ( mesh->nt ) {
1079 if(!bin) {
1080 strcpy(&chaine[0],"\n\nTriangles\n");
1081 fprintf(inm,"%s",chaine);
1082 fprintf(inm,"%" MMG5_PRId "\n",nt);
1083 } else {
1084 binch = 6; //Triangles
1085 fwrite(&binch,MMG5_SW,1,inm);
1086 bpos += (3+4*nt)*MMG5_SW; //Pos
1087 fwrite(&bpos,MMG5_SW,1,inm);
1088 fwrite(&nt,MMG5_SW,1,inm);
1089 }
1090
1091 ntreq = 0;
1092 for (k=1; k<=mesh->nt; k++) {
1093 pt = &mesh->tria[k];
1094 if ( MG_EOK(pt) ) {
1095 if(!bin) {
1096 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",mesh->point[pt->v[0]].tmp,mesh->point[pt->v[1]].tmp
1097 ,mesh->point[pt->v[2]].tmp,MMG5_abs(pt->ref));
1098 } else {
1099 fwrite(&mesh->point[pt->v[0]].tmp,MMG5_SW,1,inm);
1100 fwrite(&mesh->point[pt->v[1]].tmp,MMG5_SW,1,inm);
1101 fwrite(&mesh->point[pt->v[2]].tmp,MMG5_SW,1,inm);
1102 pt->ref = MMG5_abs(pt->ref);
1103 fwrite(&pt->ref,MMG5_SW,1,inm);
1104 }
1105
1106 if ( (pt->tag[0] & MG_REQ) && (pt->tag[1] & MG_REQ) && (pt->tag[2] & MG_REQ) ) {
1107 ntreq++;
1108 }
1109 }
1110 }
1111 if ( ntreq ) {
1112 if(!bin) {
1113 strcpy(&chaine[0],"\n\nRequiredTriangles\n");
1114 fprintf(inm,"%s",chaine);
1115 fprintf(inm,"%" MMG5_PRId "\n",ntreq);
1116 } else {
1117 binch = 17; //ReqTriangles
1118 fwrite(&binch,MMG5_SW,1,inm);
1119 bpos += (3+ntreq)*MMG5_SW; //Pos
1120 fwrite(&bpos,MMG5_SW,1,inm);
1121 fwrite(&ntreq,MMG5_SW,1,inm);
1122 }
1123 for (k=0; k<=mesh->nt; k++) {
1124 pt = &mesh->tria[k];
1125 if ( (pt->tag[0] & MG_REQ) && (pt->tag[1] & MG_REQ) && pt->tag[2] & MG_REQ ) {
1126 if(!bin) {
1127 fprintf(inm,"%" MMG5_PRId "\n",k);
1128 } else {
1129 fwrite(&k,MMG5_SW,1,inm);
1130 }
1131 }
1132 }
1133 }
1134 }
1135
1136 if ( (!mesh->xp) || (!mesh->xpoint) ) nn = ng = 0;
1137
1138 /* write normals */
1139 if ( nn ) {
1140 if(!bin) {
1141 strcpy(&chaine[0],"\n\nNormals\n");
1142 fprintf(inm,"%s",chaine);
1143 fprintf(inm,"%" MMG5_PRId "\n",nn);
1144 } else {
1145 binch = 60; //normals
1146 fwrite(&binch,MMG5_SW,1,inm);
1147 bpos += (3+(3*mesh->ver)*nn)*MMG5_SW; //Pos
1148 fwrite(&bpos,MMG5_SW,1,inm);
1149 fwrite(&nn,MMG5_SW,1,inm);
1150 }
1151 for (k=1; k<=mesh->np; k++) {
1152 ppt = &mesh->point[k];
1153 if ( !MG_VOK(ppt) ) continue;
1154 else if ( !(ppt->tag & MG_GEO) && !(ppt->tag & MG_CRN) ) {
1155 if ( ppt->tag & MG_REF ) {
1156 go = &mesh->xpoint[ppt->xp];
1157 if(!bin) {
1158 fprintf(inm,"%.15lg %.15lg %.15lg \n",go->n1[0],go->n1[1],go->n1[2]);
1159 } else {
1160 fwrite((unsigned char*)&go->n1[0],MMG5_SD,1,inm);
1161 fwrite((unsigned char*)&go->n1[1],MMG5_SD,1,inm);
1162 fwrite((unsigned char*)&go->n1[2],MMG5_SD,1,inm);
1163 }
1164 }
1165 else {
1166 if(!bin) {
1167 fprintf(inm,"%.15lg %.15lg %.15lg \n",ppt->n[0],ppt->n[1],ppt->n[2]);
1168 } else {
1169 fwrite((unsigned char*)&ppt->n[0],MMG5_SD,1,inm);
1170 fwrite((unsigned char*)&ppt->n[1],MMG5_SD,1,inm);
1171 fwrite((unsigned char*)&ppt->n[2],MMG5_SD,1,inm);
1172 }
1173 }
1174 }
1175 }
1176 if(!bin) {
1177 strcpy(&chaine[0],"\n\nNormalAtVertices\n");
1178 fprintf(inm,"%s",chaine);
1179 fprintf(inm,"%" MMG5_PRId "\n",nn);
1180 } else {
1181 binch = 20; //normalatvertices
1182 fwrite(&binch,MMG5_SW,1,inm);
1183 bpos += (3 + 2*nn)*MMG5_SW;//Pos
1184 fwrite(&bpos,MMG5_SW,1,inm);
1185 fwrite(&nn,MMG5_SW,1,inm);
1186 }
1187 nn = 0;
1188 for (k=1; k<=mesh->np; k++) {
1189 ppt = &mesh->point[k];
1190 if ( MG_VOK(ppt) && !((ppt->tag & MG_GEO) || (ppt->tag & MG_CRN) ) ) {
1191 if(!bin) {
1192 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId "\n",ppt->tmp,++nn);
1193 } else {
1194 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1195 ++nn;
1196 fwrite(&nn,MMG5_SW,1,inm);
1197 }
1198 }
1199 }
1200/*
1201 nn = 0;
1202 for (k=1; k<=mesh->nt; k++) {
1203 pt = &mesh->tria[k];
1204 if ( !MG_EOK(pt) ) continue;
1205 for (i=0; i<3; i++) {
1206 ppt = &mesh->point[pt->v[i]];
1207 if ( ppt->tag & MG_GEO ) nn++;
1208 }
1209 }
1210 GmfSetKwd(outm,GmfNormalAtTriangleVertices,nn);
1211 for (k=1; k<=mesh->nt; k++) {
1212 pt = &mesh->tria[k];
1213 if ( !MG_EOK(pt) ) continue;
1214 for (i=0; i<3; i++) {
1215 ppt = &mesh->point[pt->v[i]];
1216 if ( ppt->tag & MG_GEO ) nn++;
1217 }
1218*/
1219 }
1220
1221 /* write tangents (only if we have already analyzed the mesh =>
1222 * mesh->xpoint is allocated ) */
1223 if ( ng && mesh->xpoint ) {
1224 /* Write tangents */
1225 if(!bin) {
1226 strcpy(&chaine[0],"\n\nTangents\n");
1227 fprintf(inm,"%s",chaine);
1228 fprintf(inm,"%" MMG5_PRId "\n",ng);
1229 } else {
1230 binch = 59; //tangent
1231 fwrite(&binch,MMG5_SW,1,inm);
1232 bpos += (3+(3*mesh->ver)*ng)*MMG5_SW; //Pos
1233 fwrite(&bpos,MMG5_SW,1,inm);
1234 fwrite(&ng,MMG5_SW,1,inm);
1235 }
1236 for (k=1; k<=mesh->np; k++) {
1237 ppt = &mesh->point[k];
1238 if ( MG_VOK(ppt) && MG_EDG(ppt->tag) ) {
1239 if(!bin) {
1240 fprintf(inm,"%.15lg %.15lg %.15lg \n",ppt->n[0],ppt->n[1],ppt->n[2]);
1241 } else {
1242 fwrite((unsigned char*)&ppt->n[0],MMG5_SD,1,inm);
1243 fwrite((unsigned char*)&ppt->n[1],MMG5_SD,1,inm);
1244 fwrite((unsigned char*)&ppt->n[2],MMG5_SD,1,inm);
1245 }
1246 }
1247 }
1248 if(!bin) {
1249 strcpy(&chaine[0],"\n\nTangentAtVertices\n");
1250 fprintf(inm,"%s",chaine);
1251 fprintf(inm,"%" MMG5_PRId "\n",ng);
1252 } else {
1253 binch = 61; //tangentatvertices
1254 fwrite(&binch,MMG5_SW,1,inm);
1255 bpos += (3 + 2*ng)*MMG5_SW;//Pos
1256 fwrite(&bpos,MMG5_SW,1,inm);
1257 fwrite(&ng,MMG5_SW,1,inm);
1258 }
1259 ng = 0;
1260 for (k=1; k<=mesh->np; k++) {
1261 ppt = &mesh->point[k];
1262 if ( MG_VOK(ppt) && MG_EDG(ppt->tag) ) {
1263 if(!bin) {
1264 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId "\n",ppt->tmp,++ng);
1265 } else {
1266 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1267 ++ng;
1268 fwrite(&(ng),MMG5_SW,1,inm);
1269 }
1270 }
1271 }
1272 }
1273
1274 if ( abs(mesh->info.imprim) > 4 ) {
1275 fprintf(stdout," NUMBER OF VERTICES %8" MMG5_PRId " CORNERS %6" MMG5_PRId "\n",np,nc);
1276 fprintf(stdout," NUMBER OF TRIANGLES %8" MMG5_PRId "\n",nt);
1277
1278 if ( mesh->na )
1279 fprintf(stdout," NUMBER OF EDGES %8" MMG5_PRId " RIDGES %6" MMG5_PRId "\n",mesh->na,nr);
1280 if ( nn+ng )
1281 fprintf(stdout," NUMBER OF NORMALS %8" MMG5_PRId " TANGENTS %6" MMG5_PRId "\n",nn,ng);
1282 }
1283 /*fin fichier*/
1284 if(!bin) {
1285 strcpy(&chaine[0],"\n\nEnd\n");
1286 fprintf(inm,"%s",chaine);
1287 } else {
1288 binch = 54; //End
1289 fwrite(&binch,MMG5_SW,1,inm);
1290 bpos += 2*MMG5_SW; //bpos + End key
1291 fwrite(&bpos,MMG5_SW,1,inm);
1292 }
1293 fclose(inm);
1294 return 1;
1295}
1296
1298
1299 return MMG5_saveMshMesh(mesh,&sol,filename,1);
1300}
1301
1303
1305}
1306
1308
1309 FILE *inm;
1310 long posnp;
1311 int iswp,ier,*type,ver,bin,nsols,dim;
1312 MMG5_int k,np;
1313
1315 ier = MMG5_loadSolHeader(filename,3,&inm,&ver,&bin,&iswp,&np,&dim,&nsols,
1316 &type,&posnp,mesh->info.imprim);
1317
1318 if ( ier < 1 ) return ier;
1319
1320 if ( nsols!=1 ) {
1321 fprintf(stderr,"Error: SEVERAL SOLUTIONS FOUND (%d)\n",nsols);
1322 fclose(inm);
1323 MMG5_SAFE_FREE(type);
1324 return -1;
1325 }
1326
1327 if ( mesh->np != np ) {
1328 fprintf(stderr," ** MISMATCHES DATA: THE NUMBER OF VERTICES IN "
1329 "THE MESH (%" MMG5_PRId ") DIFFERS FROM THE NUMBER OF VERTICES IN "
1330 "THE SOLUTION (%" MMG5_PRId ") \n",mesh->np,np);
1331 fclose(inm);
1332 MMG5_SAFE_FREE(type);
1333 return -1;
1334 }
1335
1336 /* #MMG5_loadSolHeader function reads only solution at vertices so we don't
1337 have to check the entites on which the metric applies */
1338 int entities = MMG5_Vertex;
1339 ier = MMG5_chkMetricType(mesh,type,&entities,inm);
1340 if ( ier < 1 ) {
1341 MMG5_SAFE_FREE(type);
1342 return ier;
1343 }
1344
1345 /* Allocate and store the header informations for each solution */
1346 if ( !MMGS_Set_solSize(mesh,met,MMG5_Vertex,mesh->np,type[0]) ) {
1347 fclose(inm);
1348 MMG5_SAFE_FREE(type);
1349 return -1;
1350 }
1351 /* For binary file, we read the verson inside the file */
1352 if ( ver ) met->ver = ver;
1353
1354 MMG5_SAFE_FREE(type);
1355
1356 /* Read mesh solutions */
1357 rewind(inm);
1358 fseek(inm,posnp,SEEK_SET);
1359
1360 /* isotropic metric */
1361 if ( met->ver == 1 ) {
1362 /* Simple precision */
1363 for (k=1; k<=mesh->np; k++) {
1364 if ( MMG5_readFloatSol3D(met,inm,bin,iswp,k) < 0 ) return -1;
1365 }
1366 }
1367 else {
1368 /* Double precision */
1369 for (k=1; k<=mesh->np; k++) {
1370 if ( MMG5_readDoubleSol3D(met,inm,bin,iswp,k) < 0 ) return -1;
1371 }
1372 }
1373
1374 fclose(inm);
1375
1376 /* stats */
1378
1379 return 1;
1380}
1381
1383 MMG5_pSol psl;
1384 FILE *inm;
1385 long posnp;
1386 int iswp,ier,*type,ver,bin,nsols,dim;
1387 MMG5_int j,k,np;
1388 char data[16];
1389 static char mmgWarn = 0;
1390
1392 ier = MMG5_loadSolHeader(filename,3,&inm,&ver,&bin,&iswp,&np,&dim,&nsols,
1393 &type,&posnp,mesh->info.imprim);
1394 if ( ier < 1 ) return ier;
1395
1396 if ( mesh->np != np ) {
1397 fprintf(stderr," ** MISMATCHES DATA: THE NUMBER OF VERTICES IN "
1398 "THE MESH (%" MMG5_PRId ") DIFFERS FROM THE NUMBER OF VERTICES IN "
1399 "THE SOLUTION (%" MMG5_PRId ") \n",mesh->np,np);
1400 fclose(inm);
1401 MMG5_SAFE_FREE(type);
1402 return -1;
1403 }
1404
1406 mesh->nsols = nsols;
1407
1408 if ( nsols > MMG5_NSOLS_MAX ) {
1409 fprintf(stderr,"\n ## Error: %s: unexpected number of data (%d).\n",
1410 __func__,nsols);
1411 MMG5_SAFE_FREE(type);
1412 fclose(inm);
1413 return -1;
1414 }
1415
1416 if ( *sol ) MMG5_DEL_MEM(mesh,*sol);
1417 MMG5_ADD_MEM(mesh,nsols*sizeof(MMG5_Sol),"solutions array",
1418 printf(" Exit program.\n"); fclose(inm);
1419 MMG5_SAFE_FREE(type);
1420 return -1);
1421 MMG5_SAFE_CALLOC(*sol,nsols,MMG5_Sol,return -1);
1422
1423 for ( j=0; j<nsols; ++j ) {
1424 psl = *sol + j;
1425
1426 /* Give an arbitrary name to the solution because the Medit format has non
1427 * name field */
1428 sprintf(data,"sol_%" MMG5_PRId "",j);
1429 if ( !MMGS_Set_inputSolName(mesh,psl,data) ) {
1430 if ( !mmgWarn ) {
1431 mmgWarn = 1;
1432 fprintf(stderr,"\n ## Warning: %s: unable to set solution name for"
1433 " at least 1 solution.\n",__func__);
1434 }
1435 }
1436
1437 /* Allocate and store the header informations for each solution */
1438 if ( !MMGS_Set_solSize(mesh,psl,MMG5_Vertex,mesh->np,type[j]) ) {
1439 MMG5_SAFE_FREE(type);
1440 fclose(inm);
1441 return -1;
1442 }
1443 /* For binary file, we read the verson inside the file */
1444 if ( ver ) psl->ver = ver;
1445 }
1446 MMG5_SAFE_FREE(type);
1447
1448 /* read mesh solutions */
1449 rewind(inm);
1450 fseek(inm,posnp,SEEK_SET);
1451
1452 if ( (*sol)[0].ver == 1 ) {
1453 /* Simple precision */
1454 for (k=1; k<=mesh->np; k++) {
1455 for ( j=0; j<nsols; ++j ) {
1456 psl = *sol + j;
1457 if ( MMG5_readFloatSol3D(psl,inm,bin,iswp,k) < 0 ) return -1;
1458 }
1459 }
1460 }
1461 else {
1462 /* Double precision */
1463 for (k=1; k<=mesh->np; k++) {
1464 for ( j=0; j<nsols; ++j ) {
1465 psl = *sol + j;
1466 if ( MMG5_readDoubleSol3D(psl,inm,bin,iswp,k) < 0 ) return -1;
1467 }
1468 }
1469 }
1470 fclose(inm);
1471
1472 /* stats */
1474
1475 return 1;
1476}
1477
1479 FILE* inm;
1480 MMG5_pPoint ppt;
1481 int binch,bin,ier;
1482 MMG5_int bpos,k;
1483
1484 if ( !met->m ) {
1485 fprintf(stderr,"\n ## Warning: %s: no metric data to save.\n",__func__);
1486 return 1;
1487 }
1488
1489 met->ver = 2;
1490 bpos = 0;
1491 ier = MMG5_saveSolHeader( mesh,filename,&inm,met->ver,&bin,&bpos,mesh->np,met->dim,
1492 1,&met->entities,&met->type,&met->size);
1493
1494 if ( ier < 1 ) return ier;
1495
1496 for (k=1; k<=mesh->np; k++) {
1497 ppt = &mesh->point[k];
1498 if ( !MG_VOK(ppt) ) continue;
1499
1500 MMG5_writeDoubleSol3D(mesh,met,inm,bin,k,1);
1501 fprintf(inm,"\n");
1502 }
1503
1504 /*fin fichier*/
1505 if(!bin) {
1506 fprintf(inm,"\n\nEnd\n");
1507 } else {
1508 binch = 54; //End
1509 fwrite(&binch,MMG5_SW,1,inm);
1510 }
1511 fclose(inm);
1512 return 1;
1513}
1514
1516 MMG5_pSol psl;
1517 FILE* inm;
1518 MMG5_pPoint ppt;
1519 int binch,bin,ier,npointSols,ncellSols;
1520 int *type,*entities,j,*size;
1521 MMG5_int k,bpos;
1522
1523 if ( !(*sol)[0].m ) return -1;
1524
1525 (*sol)[0].ver = 2;
1526
1527 MMG5_SAFE_CALLOC(type,mesh->nsols,int,return 0);
1528 MMG5_SAFE_CALLOC(size,mesh->nsols,int,MMG5_SAFE_FREE(type);return 0);
1529 MMG5_SAFE_CALLOC(entities,mesh->nsols,int,
1530 MMG5_SAFE_FREE(type);MMG5_SAFE_FREE(size);return 0);
1531
1532 npointSols = 0;
1533 ncellSols = 0;
1534 for (k=0; k<mesh->nsols; ++k ) {
1535
1536 if ( ((*sol)[k].entities==MMG5_Noentity) || ((*sol)[k].entities==MMG5_Vertex) ) {
1537 ++npointSols;
1538 }
1539 else if ( (*sol)[k].entities == MMG5_Triangle ) {
1540 ++ncellSols;
1541 }
1542 else {
1543 printf("\n ## Warning: %s: unexpected entity type for solution %" MMG5_PRId ": %s."
1544 "\n Ignored.\n",
1545 __func__,k,MMG5_Get_entitiesName((*sol)[k].entities));
1546 }
1547
1548 type[k] = (*sol)[k].type;
1549 size[k] = (*sol)[k].size;
1550 entities[k] = (*sol)[k].entities;
1551 }
1552
1553 bpos = 0;
1554 ier = MMG5_saveSolHeader( mesh,filename,&inm,(*sol)[0].ver,&bin,&bpos,mesh->np,
1555 (*sol)[0].dim,mesh->nsols,entities,type,size);
1556
1557 if ( ier < 1 ) return ier;
1558
1559 for (k=1; k<=mesh->np; k++) {
1560 ppt = &mesh->point[k];
1561 if ( !MG_VOK(ppt) ) continue;
1562
1563 for ( j=0; j<mesh->nsols; ++j ) {
1564 psl = *sol + j;
1565
1566 if ( (psl->entities==MMG5_Noentity) || (psl->entities==MMG5_Vertex) ) {
1567 MMG5_writeDoubleSol3D(mesh,psl,inm,bin,k,0);
1568 }
1569 }
1570 fprintf(inm,"\n");
1571 }
1572
1573 MMG5_saveSolAtTrianglesHeader( mesh,inm,(*sol)[0].ver,bin,&bpos,mesh->nsols,
1574 ncellSols,entities,type,size );
1575
1576 for (k=1; k<=mesh->nt; k++) {
1577 MMG5_pTria ptt = &mesh->tria[k];
1578 if ( !MG_EOK(ptt) ) continue;
1579
1580 for ( j=0; j<mesh->nsols; ++j ) {
1581 psl = *sol + j;
1582 if ( psl->entities==MMG5_Triangle ) {
1583 MMG5_writeDoubleSol3D(mesh,psl,inm,bin,k,0);
1584 }
1585 }
1586 fprintf(inm,"\n");
1587 }
1588
1589 MMG5_SAFE_FREE(type);
1590 MMG5_SAFE_FREE(size);
1591 MMG5_SAFE_FREE(entities);
1592
1593 /* End file */
1594 if(!bin) {
1595 fprintf(inm,"\n\nEnd\n");
1596 } else {
1597 binch = 54; //End
1598 fwrite(&binch,MMG5_SW,1,inm);
1599 }
1600 fclose(inm);
1601 return 1;
1602}
1603
1605 int ier=0;
1606 const char *filenameptr,*solnameptr;
1607 char *tmp,*soltmp;
1608
1609 if ( filename && strlen(filename) ) {
1610 filenameptr = filename;
1611 solnameptr = filename;
1612 }
1613 else if (mesh->namein && strlen(mesh->namein) ) {
1614 filenameptr = mesh->namein;
1615 if ( sol && strlen(sol->namein) ) {
1616 solnameptr = sol->namein;
1617 }
1618 else {
1619 solnameptr = mesh->namein;
1620 }
1621 }
1622 else {
1623 fprintf(stderr," ## Error: %s: please provide input file name"
1624 " (either in the mesh structure or as function argument).\n",
1625 __func__);
1626 return 0;
1627 }
1628
1629 MMG5_SAFE_MALLOC(tmp,strlen(filenameptr)+1,char,return 0);
1630 strcpy(tmp,filenameptr);
1631
1632 /* read mesh/sol files */
1633 char *ptr = MMG5_Get_filenameExt(tmp);
1634 int fmt = MMG5_Get_format(ptr,MMG5_FMT_MeditASCII);
1635
1636 int8_t savesolFile = 0;
1637
1638 switch ( fmt ) {
1639 case ( MMG5_FMT_GmshASCII ): case ( MMG5_FMT_GmshBinary ):
1641 break;
1642 case ( MMG5_FMT_VtkVtp ):
1644 break;
1645 case ( MMG5_FMT_VtkVtu ):
1647 break;
1648 case ( MMG5_FMT_VtkVtk ):
1650 break;
1651 default:
1653 savesolFile = 1;
1654 break;
1655 }
1656
1657 if ( ier && savesolFile ) {
1658 /* Medit or tetgen output: save the solution in a .sol file */
1659 if ( sol && sol->np ) {
1660 MMG5_SAFE_MALLOC(soltmp,strlen(solnameptr)+1,char,return 0);
1661 strcpy(soltmp,solnameptr);
1662
1663 if ( MMGS_saveSol(mesh,sol,soltmp) == -1) {
1664 fprintf(stderr,"\n ## ERROR: WRONG DATA TYPE OR WRONG SOLUTION NUMBER.\n");
1665 ier = 0;
1666 }
1667 MMG5_SAFE_FREE(soltmp);
1668 }
1669 }
1670
1671 return ier;
1672}
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 MMGS_Set_inputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solin)
int MMGS_Set_solSize(MMG5_pMesh mesh, MMG5_pSol sol, int typEntity, MMG5_int np, int typSol)
int ier
tmp[*strlen0]
MMG5_pMesh MMG5_pSol * sol
MMG5_pMesh * mesh
MMG5_pMesh char * filename
int MMG5_readDoubleSol3D(MMG5_pSol sol, FILE *inm, int bin, int iswp, MMG5_int pos)
Definition: inout.c:2269
int MMG5_chkMetricType(MMG5_pMesh mesh, int *type, int *entities, FILE *inm)
Definition: inout.c:2661
void MMG5_printSolStats(MMG5_pMesh mesh, MMG5_pSol *sol)
Definition: inout.c:2723
float MMG5_swapf(float sbin)
Definition: inout.c:80
int MMG5_saveSolAtTrianglesHeader(MMG5_pMesh mesh, FILE *inm, int ver, int bin, MMG5_int *bpos, int nsols, int nsolsAtTriangles, int *entities, int *type, int *size)
Definition: inout.c:2524
int MMG5_saveSolHeader(MMG5_pMesh mesh, const char *filename, FILE **inm, int ver, int *bin, MMG5_int *bpos, MMG5_int np, int dim, int nsols, int *entities, int *type, int *size)
Definition: inout.c:2387
int MMG5_loadMshMesh_part1(MMG5_pMesh mesh, const char *filename, FILE **inm, long *posNodes, long *posElts, long **posNodeData, int *bin, int *iswp, MMG5_int *nelts, int *nsols)
Definition: inout.c:278
void MMG5_writeDoubleSol3D(MMG5_pMesh mesh, MMG5_pSol sol, FILE *inm, int bin, MMG5_int pos, int metricData)
Definition: inout.c:2318
int MMG5_readFloatSol3D(MMG5_pSol sol, FILE *inm, int bin, int iswp, int pos)
Definition: inout.c:2220
double MMG5_swapd(double sbin)
Definition: inout.c:97
int MMG5_loadMshMesh_part2(MMG5_pMesh mesh, MMG5_pSol *sol, FILE **inm, const long posNodes, const long posElts, const long *posNodeData, const int bin, const int iswp, const MMG5_int nelts, const int nsols)
Definition: inout.c:668
int MMG5_saveMshMesh(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename, int metricData)
Definition: inout.c:1587
void MMG5_printMetStats(MMG5_pMesh mesh, MMG5_pSol met)
Definition: inout.c:2705
int MMG5_loadSolHeader(const char *filename, int meshDim, FILE **inm, int *ver, int *bin, int *iswp, MMG5_int *np, int *dim, int *nsols, int **type, long *posnp, int imprim)
Definition: inout.c:2080
int MMG5_swapbin(int sbin)
Definition: inout.c:42
int MMGS_loadGenericMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inout_s.c:759
int MMGS_loadSol(MMG5_pMesh mesh, MMG5_pSol met, const char *filename)
Definition: inout_s.c:1307
int MMGS_loadAllSols(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inout_s.c:1382
int MMGS_loadMshMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inout_s.c:640
int MMGS_saveMesh(MMG5_pMesh mesh, const char *filename)
Definition: inout_s.c:837
int MMGS_saveMshMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inout_s.c:1297
int MMGS_saveGenericMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inout_s.c:1604
int MMGS_saveMshMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inout_s.c:1302
int MMGS_loadMshMesh_and_allData(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inout_s.c:701
int MMGS_saveSol(MMG5_pMesh mesh, MMG5_pSol met, const char *filename)
Definition: inout_s.c:1478
int MMGS_loadMesh(MMG5_pMesh mesh, const char *filename)
Definition: inout_s.c:41
int MMGS_saveAllSols(MMG5_pMesh mesh, MMG5_pSol *sol, const char *filename)
Definition: inout_s.c:1515
int MMGS_loadVtpMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:74
int MMGS_loadVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:222
int MMGS_saveVtuMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:296
int MMGS_loadVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:148
int MMGS_saveVtkMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:326
int MMGS_saveVtpMesh(MMG5_pMesh mesh, MMG5_pSol sol, const char *filename)
Definition: inoutcpp_s.cpp:356
API headers for the mmgs library.
int MMGS_zaldy(MMG5_pMesh mesh)
Definition: zaldy_s.c:263
@ MMG5_FMT_MeditBinary
Definition: libmmgtypes.h:236
@ MMG5_FMT_VtkVtk
Definition: libmmgtypes.h:243
@ MMG5_FMT_GmshBinary
Definition: libmmgtypes.h:238
@ MMG5_FMT_GmshASCII
Definition: libmmgtypes.h:237
@ MMG5_FMT_MeditASCII
Definition: libmmgtypes.h:235
@ MMG5_FMT_VtkVtp
Definition: libmmgtypes.h:242
@ MMG5_FMT_VtkVtu
Definition: libmmgtypes.h:241
#define MMG5_NSOLS_MAX
Definition: libmmgtypes.h:181
@ MMG5_Noentity
Definition: libmmgtypes.h:223
@ MMG5_Vertex
Definition: libmmgtypes.h:224
@ MMG5_Triangle
Definition: libmmgtypes.h:226
#define MG_REQ
#define MG_GEO
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
#define MG_EOK(pt)
#define MG_NUL
#define MMG5_SD
#define MMG5_ADD_MEM(mesh, size, message, law)
#define MG_EDG(tag)
#define MMG_FREAD(ptr, size, count, stream)
#define MMG5_SW
#define MG_VOK(ppt)
#define MG_CRN
#define MMG5_SAFE_MALLOC(ptr, size, type, law)
#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)
MMG5_int b
Definition: libmmgtypes.h:306
MMG5_int ref
Definition: libmmgtypes.h:307
int16_t tag
Definition: libmmgtypes.h:310
MMG5_int a
Definition: libmmgtypes.h:306
MMG mesh structure.
Definition: libmmgtypes.h:605
MMG5_int ntmax
Definition: libmmgtypes.h:612
MMG5_int nc1
Definition: libmmgtypes.h:615
MMG5_Info info
Definition: libmmgtypes.h:651
MMG5_int ne
Definition: libmmgtypes.h:612
MMG5_pPoint point
Definition: libmmgtypes.h:641
MMG5_int npmax
Definition: libmmgtypes.h:612
MMG5_pxPoint xpoint
Definition: libmmgtypes.h:642
MMG5_int xp
Definition: libmmgtypes.h:620
MMG5_int nt
Definition: libmmgtypes.h:612
MMG5_pTria tria
Definition: libmmgtypes.h:647
MMG5_int np
Definition: libmmgtypes.h:612
MMG5_pEdge edge
Definition: libmmgtypes.h:649
MMG5_int nti
Definition: libmmgtypes.h:612
MMG5_int npi
Definition: libmmgtypes.h:612
MMG5_int nprism
Definition: libmmgtypes.h:613
MMG5_int na
Definition: libmmgtypes.h:612
char * namein
Definition: libmmgtypes.h:652
Structure to store points of a MMG mesh.
Definition: libmmgtypes.h:270
double n[3]
Definition: libmmgtypes.h:272
int16_t tag
Definition: libmmgtypes.h:284
MMG5_int tmp
Definition: libmmgtypes.h:280
double c[3]
Definition: libmmgtypes.h:271
MMG5_int xp
Definition: libmmgtypes.h:279
MMG5_int ref
Definition: libmmgtypes.h:278
int entities
Definition: libmmgtypes.h:670
double * m
Definition: libmmgtypes.h:671
int16_t tag[3]
Definition: libmmgtypes.h:342
MMG5_int ref
Definition: libmmgtypes.h:335
MMG5_int v[3]
Definition: libmmgtypes.h:334
Structure to store surface points of a MMG mesh.
Definition: libmmgtypes.h:294
double n1[3]
Definition: libmmgtypes.h:295