Mmg
Simplicial remeshers (mesh adaptation, isovalue discretization, lagrangian movement)
inout.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 "mmgcommon_private.h"
37
42int MMG5_swapbin(int sbin)
43{
44 int inv;
45 char *p_in = (char *) &sbin;
46 char *p = (char *)&inv;
47
48
49 p[0] = p_in[3];
50 p[1] = p_in[2];
51 p[2] = p_in[1];
52 p[3] = p_in[0];
53
54 return inv;
55}
56
61MMG5_int MMG5_swapbin_int(MMG5_int sbin)
62{
63 MMG5_int out;
64 char *p_in = (char *) &sbin;
65 char *p_out = (char *) &out;
66 int8_t i;
67
68 for(i=0;i<8;i++)
69 {
70 p_out[i] = p_in[7-i];
71 }
72
73 return out;
74}
75
80float MMG5_swapf(float sbin)
81{
82 float out;
83 char *p_in = (char *) &sbin;
84 char *p_out = (char *) &out;
85 p_out[0] = p_in[3];
86 p_out[1] = p_in[2];
87 p_out[2] = p_in[1];
88 p_out[3] = p_in[0];
89
90 return out;
91}
92
97double MMG5_swapd(double sbin)
98{
99 double out;
100 char *p_in = (char *) &sbin;
101 char *p_out = (char *) &out;
102 int8_t i;
103
104 for(i=0;i<8;i++)
105 {
106 p_out[i] = p_in[7-i];
107 }
108
109 return out;
110}
111
126static
127int MMG5_countBinaryElts(FILE **inm, const int nelts,const int iswp,
128 int *np, int *na, int* nt,int *nq, int *ne, int *npr)
129{
130 int typ,tagNum,i,l;
131 int k,num,idx;
132 static char mmgWarn = 0;
133
134 k = 0;
135
136 while ( k<nelts ) {
137 MMG_FREAD(&typ,MMG5_SW,1,(*inm));
138 if(iswp) typ = MMG5_swapbin(typ);
139
140 switch (typ) {
141 case 1:
142 /* Edge */
143 MMG_FREAD(&num,MMG5_SW,1,(*inm));
144 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
145 if(iswp) {
146 num = MMG5_swapbin(num);
147 tagNum = MMG5_swapbin(tagNum);
148 }
149 for ( idx=0; idx<num; ++idx ) {
150 MMG_FREAD(&i,MMG5_SW,1,(*inm));
151 for ( l=0; l<tagNum; ++l ) MMG_FREAD(&i,MMG5_SW,1,(*inm));
152 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // edge->a
153 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // edge->b
154 }
155 (*na) += num;
156 k += num;
157 break;
158 case 2:
159 /* Tria */
160 MMG_FREAD(&num,MMG5_SW,1,(*inm));
161 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
162 if(iswp) {
163 num = MMG5_swapbin(num);
164 tagNum = MMG5_swapbin(tagNum);
165 }
166 for ( idx=0; idx<num; ++idx ) {
167 MMG_FREAD(&i,MMG5_SW,1,(*inm));
168 for ( l=0; l<tagNum; ++l ) MMG_FREAD(&i,MMG5_SW,1,(*inm));
169 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // tria->v[0]
170 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // tria->v[1]
171 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // tria->v[2]
172 }
173 (*nt) += num;
174 k += num;
175 break;
176 case 3:
177 /* Quad */
178 MMG_FREAD(&num,MMG5_SW,1,(*inm));
179 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
180 if(iswp) {
181 num = MMG5_swapbin(num);
182 tagNum = MMG5_swapbin(tagNum);
183 }
184 for ( idx=0; idx<num; ++idx ) {
185 MMG_FREAD(&i,MMG5_SW,1,(*inm));
186 for ( l=0; l<tagNum; ++l ) MMG_FREAD(&i,MMG5_SW,1,(*inm));
187 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // quadra->v[0]
188 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // quadra->v[1]
189 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // quadra->v[2]
190 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // quadra->v[3]
191 }
192 (*nq) += num;
193 k += num;
194 break;
195 case 4:
196 /* Tetra */
197 MMG_FREAD(&num,MMG5_SW,1,(*inm));
198 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
199 if(iswp) {
200 num = MMG5_swapbin(num);
201 tagNum = MMG5_swapbin(tagNum);
202 }
203 for ( idx=0; idx<num; ++idx ) {
204 MMG_FREAD(&i,MMG5_SW,1,(*inm));
205 for ( l=0; l<tagNum; ++l ) MMG_FREAD(&i,MMG5_SW,1,(*inm));
206 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // tetra->v[0]
207 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // tetra->v[1]
208 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // tetra->v[2]
209 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // tetra->v[3]
210 }
211 (*ne) += num;
212 k += num;
213 break;
214 case 6:
215 /* Prism */
216 MMG_FREAD(&num,MMG5_SW,1,(*inm));
217 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
218 if(iswp) {
219 num = MMG5_swapbin(num);
220 tagNum = MMG5_swapbin(tagNum);
221 }
222 for ( idx=0; idx<num; ++idx ){
223 MMG_FREAD(&i,MMG5_SW,1,(*inm));
224 for ( l=0; l<tagNum; ++l ) MMG_FREAD(&i,MMG5_SW,1,(*inm));
225 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // prism->v[0]
226 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // prism->v[1]
227 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // prism->v[2]
228 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // prism->v[3]
229 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // prism->v[4]
230 MMG_FREAD(&i,MMG5_SW,1,(*inm)); // prism->v[5]
231 }
232 (*npr) += num;
233 k += num;
234 break;
235 case 15:
236 /* Node */
237 MMG_FREAD(&num,MMG5_SW,1,(*inm));
238 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
239 if(iswp) {
240 num = MMG5_swapbin(num);
241 tagNum = MMG5_swapbin(tagNum);
242 }
243 for ( idx=0; idx<num; ++idx ){
244 MMG_FREAD(&i,MMG5_SW,1,(*inm));
245 for ( l=0; l<tagNum; ++l ) MMG_FREAD(&i,MMG5_SW,1,(*inm));
246 MMG_FREAD(&i,MMG5_SW,1,(*inm)); //node idx
247 }
248 (*np) += num;
249 k += num;
250 break;
251 default:
252 if ( !mmgWarn ) {
253 fprintf(stderr,"\n ## Warning: %s: unexpected type of element (%d) for at"
254 " least 1 element (%d).\n",__func__,typ,k);
255 mmgWarn = 1;
256 }
257 }
258 }
259 return 1;
260}
261
279 FILE **inm,
280 long *posNodes, long *posElts,
281 long **posNodeData, int *bin, int *iswp,
282 MMG5_int *nelts,int *nsols) {
283 double dbuf[9];
284 float fbuf[9];
285 int ver,oneBin,i;
286 int k,nt,na,nq,ne,npr,np;
287 int typ,tagNum,posNodeDataSize,initPosNodeDataSize;
288 char *ptr,*data,chaine[MMG5_FILESTR_LGTH],verNum[5];
289
290 ver = oneBin = 0;
291 *posNodes = 0;
292 *posElts = 0;
293 *nelts = 0;
294 *nsols = 0;
295 *bin = 0;
296 *iswp = 0;
297 mesh->np = mesh->nt = mesh->ne = 0;
298 nt = na = nq = ne = npr = np = 0;
299
300 MMG5_SAFE_CALLOC(data,strlen(filename)+7,char,return -1);
301
302 /* Allocation of the posNodeData array: we assume that we have less than 20
303 * solutions in the file (for a greater number of sol, posNoteData is
304 * reallocated) */
305 initPosNodeDataSize = posNodeDataSize = 20;
306 MMG5_SAFE_CALLOC(*posNodeData,posNodeDataSize,long,return -1);
307
308 strcpy(data,filename);
309 ptr = strstr(data,".msh");
310 if ( !ptr ) {
311 /* data contains the filename without extension */
312 strcat(data,".mshb");
313 if (!(*inm = fopen(data,"rb")) ) {
314 ptr = strstr(data,".msh");
315 *ptr = '\0';
316 strcat(data,".msh");
317 if( !((*inm) = fopen(data,"rb")) ) {
318 MMG5_SAFE_FREE(data);
319 MMG5_SAFE_FREE(*posNodeData);
320 return 0;
321 }
322 }
323 }
324 else {
325 if( !((*inm) = fopen(data,"rb")) ) {
326 MMG5_SAFE_FREE(data);
327 MMG5_SAFE_FREE(*posNodeData);
328 return 0;
329 }
330 }
331
332 if ( mesh->info.imprim >= 0 )
333 fprintf(stdout," %%%% %s OPENED\n",data);
334 MMG5_SAFE_FREE(data);
335
336
337 /* Detection of the different fields of the file */
338 strcpy(chaine,"D");
339 while(fscanf((*inm),"%127s ",&chaine[0])!=EOF ) {
340 if(!strncmp(chaine,"$MeshFormat",strlen("$MeshFormat"))) {
341 MMG_FSCANF((*inm),"%4s %d %d ",verNum,bin,&ver);
342 mesh->ver = ver/4;
343 if ( strncmp(verNum,"2.2",3) ) {
344 fprintf(stderr,"\n ## Warning: %s: format version (%s) may be not supported."
345 " Please, use the format version 2.2.\n",__func__,verNum);
346 }
347 if ( *bin ) {
348 MMG_FREAD(&oneBin,MMG5_SW,1,(*inm));
349 if ( oneBin!=1 ) {
350 assert(oneBin==16777216);
351 *iswp=1;
352 oneBin = MMG5_swapbin(oneBin);
353 }
354 }
355 continue;
356 } else if(!strncmp(chaine,"$EndMeshFormat",strlen("EndMeshFormat"))) {
357 continue;
358 } else if(!strncmp(chaine,"$Nodes",strlen("$Nodes"))) {
359 MMG_FSCANF((*inm),"%" MMG5_PRId " ",&mesh->npi);
360 *posNodes = ftell((*inm));
361 if ( *bin ) {
362 /* Skip the binary nodes data */
363 if ( mesh->ver==1 ) {
364 for ( k=1; k<=mesh->npi; ++k ) {
365 MMG_FREAD(&i,MMG5_SW,1,(*inm));
366 MMG_FREAD( &fbuf[0],MMG5_SW,3,(*inm) );
367 }
368 }
369 else {
370 for ( k=1; k<=mesh->npi; ++k ) {
371 MMG_FREAD(&i,MMG5_SW,1,(*inm));
372 MMG_FREAD( &dbuf[0],MMG5_SD,3,(*inm) );
373 }
374 }
375 }
376 continue;
377 } else if(!strncmp(chaine,"$EndNodes",strlen("$EndNodes"))) {
378 continue;
379 } else if(!strncmp(chaine,"$NodeData",strlen("$NodeData"))) {
380
381 (*posNodeData)[*nsols] = ftell((*inm));
382 if ( ++(*nsols) == posNodeDataSize ) {
383 MMG5_SAFE_RECALLOC(*posNodeData,*nsols,
384 posNodeDataSize+initPosNodeDataSize,
385 long,"posNodeData",return -1);
386 posNodeDataSize += initPosNodeDataSize;
387 }
388
389 if ( *bin ) {
390 /* Skip the binary nodes data */
391 /* String tags */
392 MMG_FSCANF((*inm),"%d ",&tagNum);
393 for ( k=0; k<tagNum; ++k ) {
394 if ( 0 != fscanf((*inm),"%*[^\n]%*c") ) {
395 fputs ( "Reading error", stderr );
396 return -1;
397 }
398 }
399
400 /* Real tags */
401 MMG_FSCANF((*inm),"%d ",&tagNum);
402 for ( k=0; k<tagNum; ++k ) {
403 if ( 0 != fscanf((*inm),"%*[^\n]%*c") ) {
404 fputs ( "Reading error", stderr );
405 return -1;
406 }
407 }
408
409 /* Integer tags */
410 MMG_FSCANF((*inm),"%d ",&tagNum);
411 MMG_FSCANF((*inm),"%d ",&i); //time step
412 MMG_FSCANF((*inm),"%d ",&typ); //type of solution
413 MMG_FSCANF((*inm),"%d ",&np);
414 for ( k=3; k<tagNum; ++k ) {
415 MMG_FSCANF((*inm),"%d",&i);
416 }
417 if ( mesh->ver==1 ) {
418 for ( k=1; k<=np; ++k ) {
419 MMG_FREAD(&i,MMG5_SW,1,(*inm));
420 MMG_FREAD( &fbuf[0],MMG5_SW,typ,(*inm) );
421 }
422 }
423 else {
424 for ( k=1; k<=np; ++k ) {
425 MMG_FREAD(&i,MMG5_SW,1,(*inm));
426 MMG_FREAD( &dbuf[0],MMG5_SD,typ,(*inm) );
427 }
428 }
429 }
430 continue;
431 } else if(!strncmp(chaine,"$EndNodeData",strlen("$EndNodeData"))) {
432 continue;
433 } else if(!strncmp(chaine,"$Elements",strlen("$Elements"))) {
434 MMG_FSCANF((*inm),"%" MMG5_PRId " ",nelts);
435 *posElts = ftell((*inm));
436
437 /* Count the elements */
438 if ( !*bin ) {
439 for ( k=0; k<*nelts; ++k) {
440 MMG_FSCANF((*inm),"%d %d ",&i,&typ);
441 switch (typ) {
442 case 1:
443 /* Edge */
444 ++na;
445 break;
446 case 2:
447 /* Tria */
448 ++nt;
449 break;
450 case 3:
451 /* Quad */
452 ++nq;
453 break;
454 case 4:
455 /* Tetra */
456 ++ne;
457 break;
458 case 6:
459 /* Prism */
460 ++npr;
461 break;
462 case 15:
463 /* Node */
464 ++np;
465 break;
466 }
467 if ( 0 != fscanf((*inm),"%*[^\n]%*c") ) {
468 fputs ( "Reading error", stderr );
469 return -1;
470 }
471 }
472 }
473 else {
474 if ( !MMG5_countBinaryElts(inm,*nelts,*iswp,&np,&na,&nt,&nq,&ne,&npr) ) {
475 fclose(*inm);
476 MMG5_SAFE_FREE(*posNodeData);
477 return -1;
478 }
479 }
480 continue;
481 }
482 else if(!strncmp(chaine,"$EndElements",strlen("$EndElements"))) {
483 continue;
484 }
485 }
486
487 if ( !mesh->npi ) {
488 fprintf(stderr," ** MISSING DATA.\n");
489 fprintf(stderr," Check that your mesh contains points and elements.\n");
490 fprintf(stderr," Exit program.\n");
491 fclose(*inm);
492 MMG5_SAFE_FREE(*posNodeData);
493 return -1;
494 }
495
496 /* memory allocation */
497 mesh->np = mesh->npi;
498 mesh->nt = mesh->nti = nt;
499 mesh->ne = mesh->nei = ne;
500 mesh->na = mesh->nai = na;
501 mesh->nprism = npr;
502 mesh->nquad = nq;
503
504 if ( !mesh->np ) {
505 fprintf(stderr," ** MISSING DATA.\n");
506 fprintf(stderr," Check that your mesh contains points.\n");
507 fprintf(stderr," Exit program.\n");
508 fclose(*inm);
509 MMG5_SAFE_FREE(*posNodeData);
510 return -1;
511 }
512
513 return 1;
514}
515
526int MMG5_check_readedMesh ( MMG5_pMesh mesh, MMG5_int nref ) {
527 MMG5_pPoint ppt;
528 MMG5_pTria ptt;
529 MMG5_pQuad pq;
530 MMG5_pPrism pp;
531 MMG5_pTetra pt;
532 int i;
533 MMG5_int k,aux;
534
535 if ( nref ) {
536 fprintf(stdout,"\n $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n");
537 fprintf(stdout," WARNING : %" MMG5_PRId " entities with unexpected refs (ref< 0).",nref);
538 fprintf(stdout," We take their absolute values.\n");
539 fprintf(stdout," $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n\n");
540 }
541
542 /* Cross the elements, mark their vertices as used and reorient the elements
543 * with negative orientation */
544 mesh->xt = 0;
545
546 if ( mesh->dim == 2 ) {
547 for (k=1; k<=mesh->nt; k++) {
548 ptt = &mesh->tria[k];
549
550 for (i=0; i<3; i++) {
551 ppt = &mesh->point[ ptt->v[i] ];
552 ppt->tag &= ~MG_NUL;
553 }
554
555 /* Set the elements references to 0 in iso mode */
556 if ( mesh->info.iso ) ptt->ref = 0;
557
558 for(i=0 ; i<3 ; i++)
559 ptt->edg[i] = 0;
560
561 if ( MMG2D_quickarea(mesh->point[ptt->v[0]].c,mesh->point[ptt->v[1]].c,
562 mesh->point[ptt->v[2]].c) < 0.0 ) {
563 /* mesh->xt temporary used to count reoriented tetra*/
564 mesh->xt++;
565 aux = ptt->v[2];
566 ptt->v[2] = ptt->v[1];
567 ptt->v[1] = aux;
568 }
569 }
570 for (k=1; k<=mesh->nquad; k++) {
571 pq = &mesh->quadra[k];
572 for (i=0; i<4; i++) {
573 ppt = &mesh->point[ pq->v[i] ];
574 ppt->tag &= ~MG_NUL;
575 }
576 }
577 }
578 else {
579 if ( mesh->ne ) {
580 /* mmg3d */
581 for ( k=1; k<=mesh->ne; k++ ) {
582 pt = &mesh->tetra[k];
583 if ( !MG_EOK(pt) ) continue;
584
585 for (i=0; i<4; i++) {
586 ppt = &mesh->point[pt->v[i]];
587 ppt->tag &= ~MG_NUL;
588 }
589
590 /* Set the elements references to 0 in iso mode */
591 if ( mesh->info.iso ) pt->ref = 0;
592
593 /* Possibly switch 2 vertices number so that each tet is positively oriented */
594 if ( MMG5_orvol(mesh->point,pt->v) < 0.0 ) {
595 /* mesh->xt temporary used to count reoriented tetra*/
596 mesh->xt++;
597 aux = pt->v[2];
598 pt->v[2] = pt->v[3];
599 pt->v[3] = aux;
600 }
601 }
602 }
603 else {
604 /* mmgs */
605 for ( k=1; k<=mesh->nt; k++ ) {
606 ptt = &mesh->tria[k];
607 if ( !MG_EOK(ptt) ) continue;
608 for (i=0; i<3; i++) {
609 ppt = &mesh->point[ptt->v[i]];
610 ppt->tag &= ~MG_NUL;
611 }
612 }
613 }
614 }
615
616 if(mesh->xt) {
617 fprintf(stdout,"\n $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n");
618 fprintf(stdout," BAD ORIENTATION : vol < 0 -- %8" MMG5_PRId " element(s) reoriented\n",mesh->xt);
619 fprintf(stdout," $$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$$ \n\n");
620 }
621 mesh->xt = 0;
622
623 /* Cross the prisms and mark their vertices as used */
624 for ( k=1; k<=mesh->nprism; k++ ) {
625 pp = &mesh->prism[k];
626 for (i=0; i<6; i++) {
627 ppt = &mesh->point[pp->v[i]];
628 ppt->tag &= ~MG_NUL;
629 }
630 }
631
632 /* stats */
633 if ( abs(mesh->info.imprim) > 3 ) {
634 fprintf(stdout," NUMBER OF VERTICES %8" MMG5_PRId "\n",mesh->np);
635 if ( mesh->ne )
636 fprintf(stdout," NUMBER OF TETRAHEDRA %8" MMG5_PRId "\n",mesh->ne);
637
638 if ( mesh->nprism )
639 fprintf(stdout," NUMBER OF PRISMS %8" MMG5_PRId "\n",mesh->nprism);
640 if ( mesh->nt )
641 fprintf(stdout," NUMBER OF TRIANGLES %8" MMG5_PRId "\n",mesh->nt);
642 if ( mesh->nquad )
643 fprintf(stdout," NUMBER OF QUADRILATERALS %8" MMG5_PRId "\n",mesh->nquad);
644 if ( mesh->na ) {
645 fprintf(stdout," NUMBER OF EDGES %8" MMG5_PRId "\n",mesh->na);
646 }
647 }
648 return 1;
649}
650
669 const long posNodes,const long posElts,
670 const long *posNodeData,const int bin,const int iswp,
671 const MMG5_int nelts,const int nsols) {
672 MMG5_pTetra pt;
673 MMG5_pPrism pp;
674 MMG5_pTria ptt;
675 MMG5_pQuad pq1;
676 MMG5_pEdge pa;
677 MMG5_pPoint ppt;
678 MMG5_pSol psl;
679 double dbuf[9];
680 float fbuf[9],fc;
681 int i,ier;
682 int l,nref,nbl_t,nbl_a,k,v[4],nt,na,nq,ne,npr,idx,ref,num,iadr;
683 int typ,tagNum;
684 int isol;
685 int8_t metricData;
686 char chaine[MMG5_FILESTR_LGTH],*ptr;
687 static int8_t mmgWarn=0, mmgWarn1=0;
688
690 rewind((*inm));
691 fseek((*inm),posNodes,SEEK_SET);
692
693 if ( mesh->ver < 2 ) {
694 for ( k=0; k< mesh->np; ++k)
695 {
696 if ( !bin ) {
697 MMG_FSCANF((*inm),"%d ",&idx);
698 ppt = &mesh->point[idx];
699 for (i=0 ; i<mesh->dim ; i++) {
700 MMG_FSCANF((*inm),"%f ",&fc);
701 ppt->c[i] = (double) fc;
702 }
703 }
704 else {
705 MMG_FREAD(&idx,MMG5_SW,1,(*inm));
706 if ( iswp ) idx = MMG5_swapbin(idx);
707 ppt = &mesh->point[idx];
708 for (i=0 ; i<mesh->dim ; i++) {
709 MMG_FREAD(&fc,MMG5_SW,1,(*inm));
710 if(iswp) fc=MMG5_swapf(fc);
711 ppt->c[i] = (double) fc;
712 }
713 }
714 ppt->tag = MG_NUL;
715 ppt->tmp = 0;
716 ppt->ref = 0;
717 }
718 }
719 else {
720 for ( k=0; k< mesh->np; ++k)
721 {
722 if ( !bin ) {
723 MMG_FSCANF((*inm),"%d ",&i);
724 ppt = &mesh->point[i];
725 MMG_FSCANF((*inm),"%lf %lf %lf ",&ppt->c[0],&ppt->c[1],&ppt->c[2]);
726 }
727 else {
728 MMG_FREAD(&i,MMG5_SW,1,(*inm));
729 if ( iswp ) i = MMG5_swapbin(i);
730 ppt = &mesh->point[i];
731 for (i=0 ; i<3 ; i++) {
732 MMG_FREAD(&ppt->c[i],MMG5_SD,1,(*inm));
733 if(iswp) ppt->c[i]=MMG5_swapd(ppt->c[i]);
734 }
735 }
736 ppt->tag = MG_NUL;
737 ppt->tmp = 0;
738 ppt->ref = 0;
739 }
740 }
741
742 rewind((*inm));
743 fseek((*inm),posElts,SEEK_SET);
744
745 nbl_a = nbl_t = nt = na = nq = ne = npr = 0;
746 nref = 0;
747
748 if ( !bin ) {
749 for ( k=0; k<nelts; ++k)
750 {
751 MMG_FSCANF((*inm),"%d %d %d ",&i,&typ, &tagNum);
752 if ( tagNum < 2 ) {
753 fprintf(stderr,"\n ## Error: %s: elt %d (type %d): Expected at least 2 tags (%d given).\n",
754 __func__,k,typ,tagNum);
755 fclose(*inm);
756 return -1;
757 }
758 MMG_FSCANF((*inm),"%d %d ",&ref,&i);
759 for ( l=2; l<tagNum; ++l ) {
760 MMG_FSCANF((*inm),"%d ",&i);
761 }
762
763 switch (typ) {
764 case 1:
765 /* Edge */
766 /* Skip edges with mesh->info.isoref refs */
767 if ( mesh->info.iso && MMG5_abs(ref)== mesh->info.isoref ) {
768 /* Skip this edge but advance the file pointer */
769 pa = &mesh->edge[0];
770 MMG_FSCANF((*inm),"%" MMG5_PRId " %" MMG5_PRId " ",&pa->a,&pa->b);
771 ++nbl_a;
772 }
773 else {
774 pa = &mesh->edge[++na];
775 MMG_FSCANF((*inm),"%" MMG5_PRId " %" MMG5_PRId " ",&pa->a,&pa->b);
776 pa->ref = ref;
777 if ( pa->ref < 0 ) {
778 pa->ref = -pa->ref;
779 ++nref;
780 }
781 pa->tag |= MG_REF;
782 }
783 assert( na+nbl_a<=mesh->na );
784 break;
785 case 2:
786 /* Tria */
787 /* Skip triangles with mesh->info.isoref refs in 3D */
788 if ( mesh->info.iso && MMG5_abs(ref)== mesh->info.isoref && mesh->dim == 3 ) {
789 /* Skip this triangle but advance the file pointer */
790 ptt = &mesh->tria[0];
791 MMG_FSCANF((*inm),"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "",&ptt->v[0],&ptt->v[1],&ptt->v[2]);
792 ++nbl_t;
793 }
794 else {
795 ptt = &mesh->tria[++nt];
796 MMG_FSCANF((*inm),"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " ",&ptt->v[0],&ptt->v[1],&ptt->v[2]);
797 ptt->ref = ref;
798 if ( ptt->ref < 0 ) {
799 ptt->ref = -ptt->ref;
800 ++nref;
801 }
802 }
803 assert( nt+nbl_t<=mesh->nt );
804 break;
805 case 3:
806 /* Quad */
807 pq1 = &mesh->quadra[++nq];
808 MMG_FSCANF((*inm),"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " ",&pq1->v[0],&pq1->v[1],&pq1->v[2],&pq1->v[3]);
809 pq1->ref = ref;
810 if ( pq1->ref < 0 ) {
811 pq1->ref = -pq1->ref;
812 ++nref;
813 }
814 assert( nq<=mesh->nquad );
815 break;
816 case 4:
817 /* Tetra for mmg3d */
818 if ( mesh->ne ) {
819 pt = &mesh->tetra[++ne];
820 MMG_FSCANF((*inm),"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " ",&pt->v[0],&pt->v[1],&pt->v[2],&pt->v[3]);
821 pt->ref = MMG5_abs(ref);
822 } else { /*skip tetra*/
823 MMG_FSCANF((*inm),"%d %d %d %d ",&v[0],&v[1],&v[2],&v[3]);
824 }
825
826 if(ref < 0) {
827 nref++;
828 }
829
830 assert( ne<=mesh->ne );
831 break;
832 case 6:
833 /* Prism for mmg3d */
834 if ( mesh->nprism )
835 {
836 pp = &mesh->prism[++npr];
837 MMG_FSCANF((*inm),"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " ",&pp->v[0],&pp->v[1],&pp->v[2],
838 &pp->v[3],&pp->v[4],&pp->v[5]);
839 pp->ref = MMG5_abs(ref);
840 }
841 if(ref < 0) {
842 nref++;
843 }
844 assert( npr<=mesh->nprism );
845 break;
846 case 15:
847 /* Node */
848 MMG_FSCANF((*inm),"%d ",&l);
849 ppt = &mesh->point[l];
850 ppt->ref = ref;
851 if ( ppt->ref < 0 ) {
852 ppt->ref = -ppt->ref;
853 ++nref;
854 }
855 assert( l<=mesh->np );
856 break;
857 default:
858 if ( !mmgWarn ) {
859 fprintf(stderr,"\n ## Warning: %s: unexpected type for at least 1 element:"
860 " element %d, type %d\n",__func__,k,typ );
861 mmgWarn = 1;
862 }
863 }
864 }
865 }
866 else {
867 k = 0;
868
869 while ( k<nelts ) {
870 MMG_FREAD(&typ,MMG5_SW,1,(*inm));
871 if(iswp) typ = MMG5_swapbin(typ);
872
873 switch (typ) {
874 case 1:
875 /* Edge */
876 MMG_FREAD(&num,MMG5_SW,1,(*inm));
877 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
878 if(iswp) {
879 num = MMG5_swapbin(num);
880 tagNum = MMG5_swapbin(tagNum);
881 }
882 if ( tagNum < 2 ) {
883 fprintf(stderr,"\n ## Error: %s: Expected at least 2 tags per element (%d given).\n",
884 __func__,tagNum);
885 fclose(*inm);
886 return -1;
887 }
888
889 for ( idx=0; idx<num; ++idx ) {
890 MMG_FREAD(&i,MMG5_SW,1,(*inm));
891
892 MMG_FREAD(&ref,MMG5_SW,1,(*inm));
893 MMG_FREAD(&i,MMG5_SW,1,(*inm));
894
895 for ( l=2; l<tagNum; ++l ) {
896 MMG_FREAD(&i,MMG5_SW,1,(*inm));
897 }
898
899 if(iswp) ref = MMG5_swapbin(ref);
900
901 /* Skip edges with mesh->info.isoref refs */
902 if ( mesh->info.iso && MMG5_abs(ref) == mesh->info.isoref ) {
903 /* Skip this edge but advance the file pointer */
904 pa = &mesh->edge[0];
905 MMG_FREAD(&pa->a,MMG5_SW,1,(*inm));
906 MMG_FREAD(&pa->b,MMG5_SW,1,(*inm));
907 ++nbl_a;
908 }
909 else {
910 pa = &mesh->edge[++na];
911 int a,b;
912 MMG_FREAD(&a,MMG5_SW,1,(*inm));
913 MMG_FREAD(&b,MMG5_SW,1,(*inm));
914 if ( iswp ) {
915 pa->a = MMG5_swapbin(a);
916 pa->b = MMG5_swapbin(b);
917 }
918 else {
919 pa->a = a;
920 pa->b = b;
921 }
922 pa->ref = ref;
923 if ( pa->ref < 0 ) {
924 pa->ref = -pa->ref;
925 ++nref;
926 }
927 pa->tag |= MG_REF;
928 }
929 assert( na+nbl_a<=mesh->na );
930 }
931 k += num;
932
933 break;
934 case 2:
935 /* Tria */
936 MMG_FREAD(&num,MMG5_SW,1,(*inm));
937 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
938 if(iswp) {
939 num = MMG5_swapbin(num);
940 tagNum = MMG5_swapbin(tagNum);
941 }
942 if ( tagNum < 2 ) {
943 fprintf(stderr,"\n ## Error: %s: Expected at least 2 tags per element (%d given).\n",
944 __func__,tagNum);
945 fclose(*inm);
946 return -1;
947 }
948
949 for ( idx=0; idx<num; ++idx ) {
950 MMG_FREAD(&i,MMG5_SW,1,(*inm));
951
952 MMG_FREAD(&ref,MMG5_SW,1,(*inm));
953 MMG_FREAD(&i,MMG5_SW,1,(*inm));
954
955 for ( l=2; l<tagNum; ++l ) {
956 MMG_FREAD(&i,MMG5_SW,1,(*inm));
957 }
958
959 if(iswp) ref = MMG5_swapbin(ref);
960
961 /* Skip triangles with mesh->info.isoref refs in 3D */
962 if ( mesh->info.iso && MMG5_abs(ref) == mesh->info.isoref && mesh->dim == 3 ) {
963 /* Skip this triangle but advance the file pointer */
964 for ( i=0; i<3 ; ++i ) {
965 MMG_FREAD(&l,MMG5_SW,1,(*inm));
966 }
967 ++nbl_t;
968 }
969 else {
970 ptt = &mesh->tria[++nt];
971 for ( i=0; i<3 ; ++i ) {
972 MMG_FREAD(&l,MMG5_SW,1,(*inm));
973 if ( iswp ) l = MMG5_swapbin(l);
974 ptt->v[i] = l;
975 }
976 ptt->ref = ref;
977 if ( ptt->ref < 0 ) {
978 ptt->ref = -ptt->ref;
979 ++nref;
980 }
981 }
982 assert( nt+nbl_t<=mesh->nt );
983 }
984 k += num;
985 break;
986 case 3:
987 /* Quad */
988 MMG_FREAD(&num,MMG5_SW,1,(*inm));
989 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
990 if(iswp) {
991 num = MMG5_swapbin(num);
992 tagNum = MMG5_swapbin(tagNum);
993 }
994 if ( tagNum < 2 ) {
995 fprintf(stderr,"\n ## Error: %s: Expected at least 2 tags per element (%d given).\n",
996 __func__,tagNum);
997 fclose(*inm);
998 return -1;
999 }
1000
1001 for ( idx=0; idx<num; ++idx ) {
1002 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1003
1004 MMG_FREAD(&ref,MMG5_SW,1,(*inm));
1005 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1006
1007 for ( l=2; l<tagNum; ++l ) {
1008 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1009 }
1010
1011 if(iswp) ref = MMG5_swapbin(ref);
1012
1013 pq1 = &mesh->quadra[++nq];
1014 for ( i=0; i<4 ; ++i ) {
1015 MMG_FREAD(&l,MMG5_SW,1,(*inm));
1016 if ( iswp ) l = MMG5_swapbin(l);
1017 pq1->v[i] = l;
1018 }
1019 pq1->ref = ref;
1020 if ( pq1->ref < 0 ) {
1021 pq1->ref = -pq1->ref;
1022 ++nref;
1023 }
1024 assert( nq<=mesh->nquad );
1025 }
1026 k += num;
1027 break;
1028 case 4:
1029 /* Tetra */
1030 MMG_FREAD(&num,MMG5_SW,1,(*inm));
1031 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
1032 if(iswp) {
1033 num = MMG5_swapbin(num);
1034 tagNum = MMG5_swapbin(tagNum);
1035 }
1036
1037 if ( tagNum < 2 ) {
1038 fprintf(stderr,"\n ## Error: %s: Expected at least 2 tags per element (%d given).\n",
1039 __func__,tagNum);
1040 fclose(*inm);
1041 return -1;
1042 }
1043
1044 for ( idx=0; idx<num; ++idx ) {
1045 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1046
1047 MMG_FREAD(&ref,MMG5_SW,1,(*inm));
1048 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1049
1050 for ( l=2; l<tagNum; ++l ) {
1051 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1052 }
1053
1054 if(iswp) ref = MMG5_swapbin(ref);
1055
1056 if ( mesh->ne ) {
1057 pt = &mesh->tetra[++ne];
1058 for ( i=0; i<4 ; ++i ) {
1059 MMG_FREAD(&l,MMG5_SW,1,(*inm));
1060 if ( iswp ) l = MMG5_swapbin(l);
1061 pt->v[i] = l;
1062 }
1063 pt->ref = MMG5_abs(ref);
1064 assert( ne<=mesh->ne );
1065 }
1066 else
1067 for ( i=0; i<4 ; ++i )
1068 MMG_FREAD(&l,MMG5_SW,1,(*inm));
1069
1070 if(ref < 0) {
1071 nref++;
1072 }
1073 }
1074 k += num;
1075 break;
1076 case 6:
1077 /* Prism */
1078 MMG_FREAD(&num,MMG5_SW,1,(*inm));
1079 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
1080 if(iswp) {
1081 num = MMG5_swapbin(num);
1082 tagNum = MMG5_swapbin(tagNum);
1083 }
1084 if ( tagNum < 2 ) {
1085 fprintf(stderr,"\n ## Error: %s: Expected at least 2 tags per element (%d given).\n",
1086 __func__,tagNum);
1087 fclose(*inm);
1088 return -1;
1089 }
1090
1091 for ( idx=0; idx<num; ++idx ) {
1092 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1093
1094 MMG_FREAD(&ref,MMG5_SW,1,(*inm));
1095 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1096
1097 for ( l=2; l<tagNum; ++l ) {
1098 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1099 }
1100
1101 if(iswp) ref = MMG5_swapbin(ref);
1102
1103 if ( mesh->nprism ) {
1104 pp = &mesh->prism[++npr];
1105 for ( i=0; i<6 ; ++i ) {
1106 MMG_FREAD(&l,MMG5_SW,1,(*inm));
1107 if ( iswp ) l = MMG5_swapbin(l);
1108 pp->v[i] = l;
1109 }
1110 pp->ref = MMG5_abs(ref);
1111 assert( npr<=mesh->nprism );
1112 }
1113 else {
1114 for ( i=0; i<6 ; ++i )
1115 MMG_FREAD(&l,MMG5_SW,1,(*inm));
1116 }
1117
1118 if(ref < 0) {
1119 nref++;
1120 }
1121 }
1122 k += num;
1123 break;
1124 case 15:
1125 /* Node */
1126 MMG_FREAD(&num,MMG5_SW,1,(*inm));
1127 MMG_FREAD(&tagNum,MMG5_SW,1,(*inm));
1128 if(iswp) {
1129 num = MMG5_swapbin(num);
1130 tagNum = MMG5_swapbin(tagNum);
1131 }
1132 if ( tagNum < 2 ) {
1133 fprintf(stderr,"\n ## Error: %s: Expected at least 2 tags per element (%d given).\n",
1134 __func__,tagNum);
1135 fclose(*inm);
1136 return -1;
1137 }
1138
1139 for ( idx=0; idx<num; ++idx ) {
1140 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1141
1142 MMG_FREAD(&ref,MMG5_SW,1,(*inm));
1143 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1144
1145 for ( l=2; l<tagNum; ++l ) {
1146 MMG_FREAD(&i,MMG5_SW,1,(*inm));
1147 }
1148
1149 if(iswp) ref = MMG5_swapbin(ref);
1150
1151 MMG_FREAD(&l,MMG5_SW,1,(*inm));
1152 if(iswp) l = MMG5_swapbin(l);
1153 ppt = &mesh->point[l];
1154 ppt->ref = ref;
1155 if ( ppt->ref < 0 ) {
1156 ppt->ref = -ppt->ref;
1157 ++nref;
1158 }
1159 assert( l<=mesh->np );
1160 }
1161 k += num;
1162 break;
1163 default:
1164 fprintf(stderr,"\n ## Error: %s: unexpected type of element (%d)\n",
1165 __func__,typ);
1166 fclose(*inm);
1167 return -1;
1168 }
1169 }
1170 }
1171
1172 /* Check data */
1173 assert ( na + nbl_a == mesh->na );
1174 assert ( nt + nbl_t == mesh->nt );
1175
1176 /* Array reallocation if ISO refs has been skipped */
1177 if ( mesh->info.iso ) {
1178 if ( mesh->nt ) {
1179 if( !nt )
1181
1182 else if ( nt < mesh->nt ) {
1183 MMG5_ADD_MEM(mesh,(nt-mesh->nt)*sizeof(MMG5_Tria),"triangles",
1184 fprintf(stderr," Exit program.\n");
1185 fclose(*inm);
1186 return -1);
1187 MMG5_SAFE_RECALLOC(mesh->tria,mesh->nt+1,(nt+1),MMG5_Tria,"triangles",
1188 return -1);
1189 }
1190 mesh->nt = nt;
1191 }
1192 if ( mesh->na ) {
1193 if( !na )
1195 else if ( na < mesh->na ) {
1196 MMG5_ADD_MEM(mesh,(na-mesh->na)*sizeof(MMG5_Edge),"edges",
1197 fprintf(stderr," Exit program.\n");
1198 fclose(*inm);
1199 return -1);
1200 MMG5_SAFE_RECALLOC(mesh->edge,mesh->na+1,(na+1),MMG5_Edge,"edges",
1201 return -1);
1202 }
1203 mesh->na = na;
1204 }
1205 }
1206
1208 if ( ier < 1 ) return ier;
1209
1210 if ( sol && *sol ) {
1212 /* Init (*sol)[0] for the case where nsols=0 */
1213 psl = *sol;
1214 psl->ver = mesh->ver;
1215 psl->dim = mesh->dim;
1216 psl->type = 1;
1217
1218 for ( isol=0; isol < nsols; ++isol ) {
1219 assert ( posNodeData[isol] );
1220
1221 rewind((*inm));
1222 fseek((*inm),posNodeData[isol],SEEK_SET);
1223
1224 psl = *sol + isol;
1225
1226 psl->ver = mesh->ver;
1227 psl->dim = mesh->dim;
1228 psl->entities = MMG5_Vertex;
1229 psl->type = 1;
1230
1231 /* String tags: The first one stores the solution name */
1232 MMG_FSCANF((*inm),"%d ",&tagNum);
1233 if ( 1 != fscanf(*inm,"\"%127s\"\n",&chaine[0]) ) {
1234 MMG_FSCANF(*inm,"%127s\n",&chaine[0]);
1235 }
1236
1237 ptr = NULL;
1238 ptr = strstr(chaine,":metric");
1239
1240 metricData = 0;
1241 if ( ptr ) {
1242 *ptr = '\0';
1243 metricData = 1;
1244 }
1245
1246
1247 if ( !MMG5_Set_inputSolName(mesh,psl,chaine) ) {
1248 if ( !mmgWarn1 ) {
1249 mmgWarn1 = 1;
1250 fprintf(stderr,"\n ## Warning: %s: unable to set solution name for"
1251 " at least 1 solution.\n",__func__);
1252 }
1253 }
1254
1255 for ( k=1; k<tagNum; ++k ) {
1256 if ( 0 != fscanf((*inm),"%*[^\n]%*c") ) return -1;
1257 }
1258
1259 /* Real tags ignored */
1260 if ( fscanf((*inm),"%d",&tagNum) ) {
1261 for ( k=0; k<tagNum; ++k ) {
1262 MMG_FSCANF((*inm),"%f",&fc);
1263 }
1264 }
1265
1266 /* Integer tags : allow to recover the number of sols and their types */
1267 MMG_FSCANF((*inm),"%d ",&tagNum);
1268 if ( tagNum < 3 ) {
1269 fprintf(stderr," Error: %s: node data: Expected at least 3 tags (%d given).\n",
1270 __func__,tagNum);
1271 fclose(*inm);
1272 return -1;
1273 }
1274
1275 MMG_FSCANF((*inm),"%d ",&i); //time step;
1276 MMG_FSCANF((*inm),"%d ",&typ); //type of solution: 1=scalar, 3=vector, 9=tensor ;
1277 MMG_FSCANF((*inm),"%" MMG5_PRId " ",&psl->np);
1278
1279 for ( k=3; k<tagNum; ++k ) {
1280 MMG_FSCANF((*inm),"%d",&i);
1281 }
1282
1283 if ( mesh->np != psl->np ) {
1284 fprintf(stderr," ** MISMATCHES DATA: THE NUMBER OF VERTICES IN "
1285 "THE MESH (%" MMG5_PRId ") DIFFERS FROM THE NUMBER OF VERTICES IN "
1286 "THE SOLUTION (%" MMG5_PRId ") \n",mesh->np,psl->np);
1287 fclose(*inm);
1288 return -1;
1289 }
1290
1291 if ( typ == 1 ) {
1292 psl->size = 1;
1293 psl->type = 1;
1294 }
1295 else if ( typ == 3 ) {
1296 psl->size = psl->dim;
1297 psl->type = 2;
1298 }
1299 else if ( typ == 9 ) {
1300 if ( metricData ) {
1301 psl->size = (psl->dim*(psl->dim+1))/2;
1302 psl->type = 3;
1303 }
1304 else {
1305 psl->size = psl->dim*psl->dim;
1306 psl->type = 4;
1307 }
1308 }
1309 else {
1310 fprintf(stderr," ** DATA TYPE IGNORED %d \n",typ);
1311 fclose(*inm);
1312 return -1;
1313 }
1314
1315 /* mem alloc */
1316 if ( psl->m ) MMG5_DEL_MEM(mesh,psl->m);
1317 psl->npmax = mesh->npmax;
1318
1319 MMG5_ADD_MEM(mesh,(psl->size*(psl->npmax+1))*sizeof(double),"initial solution",
1320 fprintf(stderr," Exit program.\n");
1321 fclose(*inm);
1322 return -1);
1323 MMG5_SAFE_CALLOC(psl->m,psl->size*(psl->npmax+1),double,return -1);
1324
1325 /* isotropic solution */
1326 if ( psl->size == 1 ) {
1327 if ( psl->ver == 1 ) {
1328 for (k=1; k<=psl->np; k++) {
1329 if(!bin){
1330 MMG_FSCANF((*inm),"%d ",&idx);
1331 MMG_FSCANF((*inm),"%f ",&fbuf[0]);
1332 } else {
1333 MMG_FREAD(&idx,MMG5_SW,1, (*inm));
1334 if(iswp) idx = MMG5_swapbin(idx);
1335 MMG_FREAD(&fbuf[0],MMG5_SW,1,(*inm));
1336 if(iswp) fbuf[0]=MMG5_swapf(fbuf[0]);
1337 }
1338 psl->m[idx] = fbuf[0];
1339 }
1340 }
1341 else {
1342 for (k=1; k<=psl->np; k++) {
1343 if(!bin){
1344 MMG_FSCANF((*inm),"%d ",&idx);
1345 MMG_FSCANF((*inm),"%lf ",&dbuf[0]);
1346 } else {
1347 MMG_FREAD(&idx,MMG5_SW,1, (*inm));
1348 if(iswp) idx = MMG5_swapbin(idx);
1349 MMG_FREAD(&dbuf[0],MMG5_SD,1,(*inm));
1350 if(iswp) dbuf[0]=MMG5_swapd(dbuf[0]);
1351 }
1352 psl->m[idx] = dbuf[0];
1353 }
1354 }
1355 }
1356 /* vector displacement only */
1357 else if ( psl->size == psl->dim ) {
1358 if ( psl->ver == 1 ) {
1359 for (k=1; k<=psl->np; k++) {
1360 if(!bin){
1361 MMG_FSCANF((*inm),"%d ",&idx);
1362 for (i=0; i<psl->dim; i++) {
1363 MMG_FSCANF((*inm),"%f ",&fbuf[0]);
1364 psl->m[psl->dim*idx+i] = fbuf[0];
1365 }
1366 } else {
1367 MMG_FREAD(&idx,MMG5_SW,1, (*inm));
1368 if(iswp) idx = MMG5_swapbin(idx);
1369 for (i=0; i<psl->dim; i++) {
1370 MMG_FREAD(&fbuf[0],MMG5_SW,1,(*inm));
1371 if(iswp) fbuf[0]=MMG5_swapf(fbuf[0]);
1372 psl->m[psl->dim*idx+i] = fbuf[0];
1373 }
1374 }
1375 }
1376 }
1377 else {
1378 for (k=1; k<=psl->np; k++) {
1379 if(!bin){
1380 MMG_FSCANF((*inm),"%d ",&idx);
1381
1382 for (i=0; i<psl->dim; i++) {
1383 MMG_FSCANF((*inm),"%lf ",&dbuf[0]);
1384 psl->m[psl->dim*idx+i] = dbuf[0];
1385 }
1386
1387 } else {
1388 MMG_FREAD(&idx,MMG5_SW,1, (*inm));
1389 if(iswp) idx = MMG5_swapbin(idx);
1390
1391 for (i=0; i<psl->dim; i++) {
1392 MMG_FREAD(&dbuf[0],MMG5_SD,1,(*inm));
1393 if(iswp) dbuf[0]=MMG5_swapd(dbuf[0]);
1394 psl->m[psl->dim*idx+i] = dbuf[0];
1395 }
1396
1397 }
1398 }
1399 }
1400 }
1401 /* anisotropic sol */
1402 else {
1403 if ( psl->ver == 1 ) {
1404 /* Solution at simple precision */
1405 for (k=1; k<=psl->np; k++) {
1406
1407 if(!bin){
1408 MMG_FSCANF((*inm),"%d ",&idx);
1409 for(i=0 ; i<9 ; i++)
1410 MMG_FSCANF((*inm),"%f ",&fbuf[i]);
1411 } else {
1412 MMG_FREAD(&idx,MMG5_SW,1, (*inm));
1413 if(iswp) idx = MMG5_swapbin(idx);
1414 for(i=0 ; i<9 ; i++) {
1415 MMG_FREAD(&fbuf[i],MMG5_SW,1,(*inm));
1416 if(iswp) fbuf[i]=MMG5_swapf(fbuf[i]);
1417 }
1418 }
1419
1420 iadr = psl->size*idx;
1421
1422 if ( !metricData ) {
1423 if ( psl->dim ==2 ) {
1424 psl->m[iadr] = fbuf[0];
1425 psl->m[iadr+1] = fbuf[1];
1426 psl->m[iadr+2] = fbuf[3];
1427 psl->m[iadr+3] = fbuf[4];
1428 }
1429 else {
1430 for(i=0 ; i<9 ; i++) {
1431 psl->m[iadr+i] = fbuf[i];
1432 }
1433 }
1434 }
1435 else {
1436 if ( psl->dim ==2 ) {
1437 assert ( fbuf[1] == fbuf[3] );
1438
1439 psl->m[iadr] = fbuf[0];
1440 psl->m[iadr+1] = fbuf[1];
1441 psl->m[iadr+2] = fbuf[4];
1442 }
1443 else {
1444 assert ( fbuf[1]==fbuf[3] && fbuf[2]==fbuf[6] && fbuf[5]==fbuf[7] );
1445
1446 psl->m[iadr+0] = fbuf[0];
1447 psl->m[iadr+1] = fbuf[1];
1448 psl->m[iadr+2] = fbuf[2];
1449 psl->m[iadr+3] = fbuf[4];
1450 psl->m[iadr+4] = fbuf[5];
1451 psl->m[iadr+5] = fbuf[8];
1452 }
1453 }
1454 }
1455 }
1456 else {
1457 for (k=1; k<=psl->np; k++) {
1458 /* Solution at double precision */
1459 if(!bin){
1460 MMG_FSCANF((*inm),"%d ",&idx);
1461 for(i=0 ; i<9 ; i++)
1462 MMG_FSCANF((*inm),"%lf ",&dbuf[i]);
1463 } else {
1464 MMG_FREAD(&idx,MMG5_SW,1, (*inm));
1465 if(iswp) idx = MMG5_swapbin(idx);
1466 for(i=0 ; i<9 ; i++) {
1467 MMG_FREAD(&dbuf[i],MMG5_SD,1,(*inm));
1468 if(iswp) dbuf[i]=MMG5_swapd(dbuf[i]);
1469 }
1470 }
1471
1472 iadr = psl->size*idx;
1473
1474 if ( !metricData ) {
1475 if ( psl->dim ==2 ) {
1476 psl->m[iadr ] = dbuf[0];
1477 psl->m[iadr+1] = dbuf[1];
1478 psl->m[iadr+2] = dbuf[3];
1479 psl->m[iadr+3] = dbuf[4];
1480 }
1481 else {
1482 for(i=0 ; i<9 ; i++) {
1483 psl->m[iadr+i] = dbuf[i];
1484 }
1485 }
1486 }
1487 else {
1488 if ( psl->dim ==2 ) {
1489 assert ( dbuf[1] == dbuf[3] );
1490
1491 psl->m[iadr ] = dbuf[0];
1492 psl->m[iadr+1] = dbuf[1];
1493 psl->m[iadr+2] = dbuf[4];
1494 }
1495 else {
1496 assert ( dbuf[1]==dbuf[3] || dbuf[2]==dbuf[6] || dbuf[5]==dbuf[7] );
1497
1498 psl->m[iadr+0] = dbuf[0];
1499 psl->m[iadr+1] = dbuf[1];
1500 psl->m[iadr+2] = dbuf[2];
1501 psl->m[iadr+3] = dbuf[4];
1502 psl->m[iadr+4] = dbuf[5];
1503 psl->m[iadr+5] = dbuf[8];
1504 }
1505 }
1506 }
1507 }
1508 }
1509
1510 psl->npi = psl->np;
1511 }
1512 }
1513
1514 fclose((*inm));
1515
1516 return 1;
1517}
1518
1529 double dbuf[6]) {
1530 MMG5_pPoint ppt;
1531 double mtmp[3],r[3][3];
1532 int i;
1533
1534 ppt = &mesh->point[ip];
1535 if ( mesh->info.metRidTyp &&
1536 ( !(MG_SIN(ppt->tag) || (ppt->tag & MG_NOM) || (ppt->tag & MG_NOSURF))
1537 && (ppt->tag & MG_GEO) ) ) {
1538 /* Specific storage of the aniso metric at ridge */
1539 if ( mesh->xp ) {
1540 // Arbitrary, we take the metric associated to the surface ruled by n_1
1541 mtmp[0] = sol->m[sol->size*(ip)];
1542 mtmp[1] = sol->m[sol->size*(ip)+1];
1543 mtmp[2] = sol->m[sol->size*(ip)+3];
1544
1545 // Rotation matrix.
1546 r[0][0] = ppt->n[0];
1547 r[1][0] = ppt->n[1];
1548 r[2][0] = ppt->n[2];
1549 r[0][1] = mesh->xpoint[ppt->xp].n1[1]*ppt->n[2]
1550 - mesh->xpoint[ppt->xp].n1[2]*ppt->n[1];
1551 r[1][1] = mesh->xpoint[ppt->xp].n1[2]*ppt->n[0]
1552 - mesh->xpoint[ppt->xp].n1[0]*ppt->n[2];
1553 r[2][1] = mesh->xpoint[ppt->xp].n1[0]*ppt->n[1]
1554 - mesh->xpoint[ppt->xp].n1[1]*ppt->n[0];
1555 r[0][2] = mesh->xpoint[ppt->xp].n1[0];
1556 r[1][2] = mesh->xpoint[ppt->xp].n1[1];
1557 r[2][2] = mesh->xpoint[ppt->xp].n1[2];
1558
1559 // Metric in the canonic space
1560 dbuf[0] = mtmp[0]*r[0][0]*r[0][0] + mtmp[1]*r[0][1]*r[0][1] + mtmp[2]*r[0][2]*r[0][2];
1561 dbuf[1] = mtmp[0]*r[0][0]*r[1][0] + mtmp[1]*r[0][1]*r[1][1] + mtmp[2]*r[0][2]*r[1][2];
1562 dbuf[2] = mtmp[0]*r[0][0]*r[2][0] + mtmp[1]*r[0][1]*r[2][1] + mtmp[2]*r[0][2]*r[2][2];
1563 dbuf[3] = mtmp[0]*r[1][0]*r[1][0] + mtmp[1]*r[1][1]*r[1][1] + mtmp[2]*r[1][2]*r[1][2];
1564 dbuf[4] = mtmp[0]*r[1][0]*r[2][0] + mtmp[1]*r[1][1]*r[2][1] + mtmp[2]*r[1][2]*r[2][2];
1565 dbuf[5] = mtmp[0]*r[2][0]*r[2][0] + mtmp[1]*r[2][1]*r[2][1] + mtmp[2]*r[2][2]*r[2][2];
1566 }
1567 else { // Cannot recover the metric
1568 for (i=0; i<sol->size; i++) dbuf[i] = 0.;
1569 }
1570 }
1571 else {
1572 for (i=0; i<sol->size; i++) dbuf[i] = sol->m[sol->size*ip+i];
1573 }
1574}
1575
1588 int metricData) {
1589 FILE* inm;
1590 MMG5_pPoint ppt;
1591 MMG5_pTetra pt;
1592 MMG5_pPrism pp;
1593 MMG5_pTria ptt;
1594 MMG5_pQuad pq;
1595 MMG5_pEdge pa;
1596 MMG5_pSol psl;
1597 double dbuf[6];
1598 int bin,i,typ;
1599 MMG5_int header[3],nq,ne,npr,np,nt,na,k,iadr,nelts,word;
1600 int isol,nsols;
1601 char *ptr,*data;
1602 static char mmgWarn = 0;
1603
1604 bin = 0;
1605
1606 MMG5_SAFE_CALLOC(data,strlen(filename)+7,char,return 0);
1607 strcpy(data,filename);
1608
1609 ptr = strstr(data,".msh");
1610 if ( !ptr ) {
1611 /* data contains the filename without extension */
1612 strcat(data,".mshb");
1613 if (!(inm = fopen(data,"wb")) ) {
1614 ptr = strstr(data,".msh");
1615 *ptr = '\0';
1616 strcat(data,".msh");
1617 if( !(inm = fopen(data,"wb")) ) {
1618 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
1619 MMG5_SAFE_FREE(data);
1620 return 0;
1621 }
1622 }
1623 else bin=1;
1624 }
1625 else {
1626 ptr = strstr(data,".mshb");
1627 if ( ptr ) bin = 1;
1628 if( !(inm = fopen(data,"wb")) ) {
1629 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
1630 MMG5_SAFE_FREE(data);
1631 return 0;
1632 }
1633 }
1634
1635 if ( mesh->info.imprim >= 0 )
1636 fprintf(stdout," %%%% %s OPENED\n",data);
1637 MMG5_SAFE_FREE(data);
1638
1639 /* Entete fichier*/
1640 fprintf(inm,"$MeshFormat\n");
1641 fprintf(inm,"2.2 %d %d\n",bin,8);
1642 if ( bin ) {
1643 word = 1;
1644 fwrite(&word,MMG5_SW,1,inm);
1645 fprintf(inm,"\n");
1646 }
1647 fprintf(inm,"$EndMeshFormat\n");
1648
1649 /* Vertices */
1650 np = 0;
1651 for (k=1; k<=mesh->np; k++) {
1652 ppt = &mesh->point[k];
1653 if ( MG_VOK(ppt) ) {
1654 ppt->tmp = ++np;
1655 if ( mesh->dim==2 ) ppt->c[2] = 0.;
1656 }
1657 }
1658 fprintf(inm,"$Nodes\n");
1659 fprintf(inm,"%" MMG5_PRId "\n",np);
1660
1661 for (k=1; k<=mesh->np; k++) {
1662 ppt = &mesh->point[k];
1663 if ( MG_VOK(ppt) ) {
1664 if(!bin) {
1665 fprintf(inm," %" MMG5_PRId "",ppt->tmp);
1666 for ( i=0; i<3; ++i )
1667 fprintf(inm," %.15lg",ppt->c[i]);
1668 fprintf(inm,"\n");
1669 } else {
1670 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1671 fwrite(&ppt->c[0],MMG5_SD,3,inm);
1672 }
1673 }
1674 }
1675 if ( bin ) fprintf(inm,"\n");
1676 fprintf(inm,"$EndNodes\n");
1677
1678 /* Data reading */
1680 ne = 0;
1681 for (k=1; k<=mesh->ne; k++) {
1682 pt = &mesh->tetra[k];
1683 if ( !MG_EOK(pt) ) continue;
1684 ++ne;
1685 }
1686 npr = 0;
1687 for (k=1; k<=mesh->nprism; k++) {
1688 pp = &mesh->prism[k];
1689 if ( !MG_EOK(pp) ) continue;
1690 npr++;
1691 }
1692 nq = 0;
1693 for (k=1; k<=mesh->nquad; k++) {
1694 pq = &mesh->quadra[k];
1695 if ( !MG_EOK(pq) ) continue;
1696 nq++;
1697 }
1698 nt = 0;
1699 for (k=1; k<=mesh->nt; k++) {
1700 ptt = &mesh->tria[k];
1701 if ( !MG_EOK(ptt) ) continue;
1702 nt++;
1703 }
1704 na = 0;
1705 for (k=1; k<=mesh->na; k++) {
1706 pa = &mesh->edge[k];
1707 if ( (!pa || !pa->a) ) continue;
1708 na++;
1709 }
1710
1711 fprintf(inm,"$Elements\n");
1712 fprintf(inm,"%" MMG5_PRId "\n", np + ne + npr + nt + nq + na );
1713
1716 nelts = 0;
1717
1718 /* Nodes */
1719 if ( bin ) {
1720 header[0] = 15;// Node keyword
1721 header[1] = np;
1722 header[2] = 2; // 2 tags per node
1723 fwrite(header, MMG5_SW, 3, inm);
1724 }
1725
1726 for ( k=1; k<= mesh->np; ++k)
1727 {
1728 ppt = &mesh->point[k];
1729 if ( !MG_VOK(ppt) ) continue;
1730 ++nelts;
1731
1732 if ( !bin ) fprintf(inm,"%" MMG5_PRId " 15 2 %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",nelts,MMG5_abs(ppt->ref),
1733 MMG5_abs(ppt->ref),ppt->tmp);
1734 else {
1735 fwrite(&nelts,MMG5_SW,1,inm);
1736 word = MMG5_abs(ppt->ref);
1737 fwrite(&word,MMG5_SW,1,inm);
1738 fwrite(&word,MMG5_SW,1,inm);
1739 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1740 }
1741 }
1742
1743 /* Edges */
1744 if ( bin && na ) {
1745 header[0] = 1;// Edge keyword
1746 header[1] = na;
1747 header[2] = 2; // 2 tags per edge
1748 fwrite(header, MMG5_SW, 3, inm);
1749 }
1750
1751 for (k=1; k<=mesh->na; ++k) {
1752 pa = &mesh->edge[k];
1753
1754 if ( !pa || !pa->a ) continue;
1755 ++nelts;
1756
1757 if(!bin)
1758 fprintf(inm,"%" MMG5_PRId " 1 2 %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",nelts,pa->ref,pa->ref,
1759 mesh->point[pa->a].tmp,mesh->point[pa->b].tmp);
1760 else {
1761 fwrite(&nelts,MMG5_SW,1,inm);
1762 fwrite(&pa->ref,MMG5_SW,1,inm);
1763 fwrite(&pa->ref,MMG5_SW,1,inm);
1764 fwrite(&mesh->point[pa->a].tmp,MMG5_SW,1,inm);
1765 fwrite(&mesh->point[pa->b].tmp,MMG5_SW,1,inm);
1766 }
1767 }
1768
1769 /* Triangles */
1770 if ( bin && nt ) {
1771 header[0] = 2;// Tria keyword
1772 header[1] = nt;
1773 header[2] = 2; // 2 tags per tria
1774 fwrite(header, MMG5_SW, 3, inm);
1775 }
1776
1777 for (k=1; k<=mesh->nt; ++k) {
1778 ptt = &mesh->tria[k];
1779
1780 if ( !MG_EOK(ptt) ) continue;
1781 ++nelts;
1782
1783 if(!bin)
1784 fprintf(inm,"%" MMG5_PRId " 2 2 %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",nelts,ptt->ref,ptt->ref,
1785 mesh->point[ptt->v[0]].tmp,mesh->point[ptt->v[1]].tmp,
1786 mesh->point[ptt->v[2]].tmp);
1787 else {
1788 fwrite(&nelts,MMG5_SW,1,inm);
1789 fwrite(&ptt->ref,MMG5_SW,1,inm);
1790 fwrite(&ptt->ref,MMG5_SW,1,inm);
1791 fwrite(&mesh->point[ptt->v[0]].tmp,MMG5_SW,1,inm);
1792 fwrite(&mesh->point[ptt->v[1]].tmp,MMG5_SW,1,inm);
1793 fwrite(&mesh->point[ptt->v[2]].tmp,MMG5_SW,1,inm);
1794 }
1795 }
1796
1797 /* Quads */
1798 if ( bin && nq ) {
1799 header[0] = 3;// Quad keyword
1800 header[1] = nq;
1801 header[2] = 2; // 2 tags per quad
1802 fwrite(header, MMG5_SW, 3, inm);
1803 }
1804
1805 for (k=1; k<=mesh->nquad; ++k) {
1806 pq = &mesh->quadra[k];
1807 if ( !MG_EOK(pq) ) continue;
1808 ++nelts;
1809
1810 if(!bin)
1811 fprintf(inm,"%" MMG5_PRId " 3 2 %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",nelts,pq->ref,pq->ref,
1812 mesh->point[pq->v[0]].tmp,mesh->point[pq->v[1]].tmp,
1813 mesh->point[pq->v[2]].tmp,mesh->point[pq->v[3]].tmp);
1814 else {
1815 fwrite(&nelts,MMG5_SW,1,inm);
1816 fwrite(&pq->ref,MMG5_SW,1,inm);
1817 fwrite(&pq->ref,MMG5_SW,1,inm);
1818 fwrite(&mesh->point[pq->v[0]].tmp,MMG5_SW,1,inm);
1819 fwrite(&mesh->point[pq->v[1]].tmp,MMG5_SW,1,inm);
1820 fwrite(&mesh->point[pq->v[2]].tmp,MMG5_SW,1,inm);
1821 fwrite(&mesh->point[pq->v[3]].tmp,MMG5_SW,1,inm);
1822 }
1823 }
1824
1825 /* Tetra */
1826 if ( bin && ne ) {
1827 header[0] = 4;// Tetra keyword
1828 header[1] = ne;
1829 header[2] = 2; // 2 tags per quad
1830 fwrite(header, MMG5_SW, 3, inm);
1831 }
1832
1833 for (k=1; k<=mesh->ne; ++k) {
1834 pt = &mesh->tetra[k];
1835 if ( !MG_EOK(pt) ) continue;
1836 ++nelts;
1837
1838 if(!bin)
1839 fprintf(inm,"%" MMG5_PRId " 4 2 %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",nelts,pt->ref,pt->ref,
1840 mesh->point[pt->v[0]].tmp,mesh->point[pt->v[1]].tmp,
1841 mesh->point[pt->v[2]].tmp,mesh->point[pt->v[3]].tmp);
1842 else {
1843 fwrite(&nelts,MMG5_SW,1,inm);
1844 fwrite(&pt->ref,MMG5_SW,1,inm);
1845 fwrite(&pt->ref,MMG5_SW,1,inm);
1846 fwrite(&mesh->point[pt->v[0]].tmp,MMG5_SW,1,inm);
1847 fwrite(&mesh->point[pt->v[1]].tmp,MMG5_SW,1,inm);
1848 fwrite(&mesh->point[pt->v[2]].tmp,MMG5_SW,1,inm);
1849 fwrite(&mesh->point[pt->v[3]].tmp,MMG5_SW,1,inm);
1850 }
1851 }
1852
1853 /* Prisms */
1854 if ( bin && npr ) {
1855 header[0] = 6;// Prism keyword
1856 header[1] = npr;
1857 header[2] = 2; // 2 tags per prism
1858 fwrite(header, MMG5_SW, 3, inm);
1859 }
1860
1861 for (k=1; k<=mesh->nprism; ++k) {
1862 pp = &mesh->prism[k];
1863 if ( !MG_EOK(pp) ) continue;
1864 ++nelts;
1865
1866 if(!bin)
1867 fprintf(inm,"%" MMG5_PRId " 6 2 %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",nelts,pp->ref,pp->ref,
1868 mesh->point[pp->v[0]].tmp,mesh->point[pp->v[1]].tmp,
1869 mesh->point[pp->v[2]].tmp,mesh->point[pp->v[3]].tmp,
1870 mesh->point[pp->v[4]].tmp,mesh->point[pp->v[5]].tmp);
1871 else {
1872 fwrite(&nelts,MMG5_SW,1,inm);
1873 fwrite(&pp->ref,MMG5_SW,1,inm);
1874 fwrite(&pp->ref,MMG5_SW,1,inm);
1875 fwrite(&mesh->point[pp->v[0]].tmp,MMG5_SW,1,inm);
1876 fwrite(&mesh->point[pp->v[1]].tmp,MMG5_SW,1,inm);
1877 fwrite(&mesh->point[pp->v[2]].tmp,MMG5_SW,1,inm);
1878 fwrite(&mesh->point[pp->v[3]].tmp,MMG5_SW,1,inm);
1879 fwrite(&mesh->point[pp->v[4]].tmp,MMG5_SW,1,inm);
1880 fwrite(&mesh->point[pp->v[5]].tmp,MMG5_SW,1,inm);
1881 }
1882 }
1883 if ( bin ) fprintf(inm,"\n");
1884 fprintf(inm,"$EndElements\n");
1885
1886 /* stats */
1887 if ( abs(mesh->info.imprim) > 3 ) {
1888 fprintf(stdout," NUMBER OF VERTICES %8" MMG5_PRId "\n",np);
1889 fprintf(stdout," NUMBER OF TETRAHEDRA %8" MMG5_PRId "\n",ne);
1890 if ( npr )
1891 fprintf(stdout," NUMBER OF PRISMS %8" MMG5_PRId "\n",npr);
1892 if ( nt )
1893 fprintf(stdout," NUMBER OF TRIANGLES %8" MMG5_PRId "\n",nt);
1894 if ( nq )
1895 fprintf(stdout," NUMBER OF QUADRILATERALS %8" MMG5_PRId "\n",nq);
1896 if ( na ) {
1897 fprintf(stdout," NUMBER OF EDGES %8" MMG5_PRId "\n",na);
1898 }
1899
1900 }
1901
1903 if ( metricData==1 ) {
1904 if ( sol && *sol && sol[0]->np ) {
1905 nsols = 1;
1906 }
1907 else {
1908 /* In analysis mode (-noinsert -noswap -nomove), metric is not allocated */
1909 nsols = 0;
1910 }
1911 }
1912 else {
1913 nsols = mesh->nsols;
1914 }
1915
1916 for ( isol=0; isol<nsols; ++isol) {
1917 psl = *sol + isol;
1918
1919 if ( !psl ) {
1920 continue;
1921 }
1922
1923 if ( !psl->m ) {
1924 if ( !mmgWarn ) {
1925 mmgWarn = 1;
1926 fprintf(stderr, " ## Warning: %s: missing data for at least 1 solution."
1927 " Skipped.\n",__func__);
1928 }
1929 continue;
1930 }
1931
1932 fprintf(inm,"$NodeData\n");
1933
1934 /* One string tag saying the type of solution saved.
1935 Add the "metric" keyword if we save a metric */
1936 fprintf(inm,"1\n");
1937
1938 if ( psl->size == 1 ) {
1939 typ = 1;
1940 }
1941 else if ( psl->size == psl->dim ) {
1942 typ = 3;
1943 }
1944 else {
1945 typ = 9;
1946 }
1947
1948 if ( psl->namein ) {
1949 char *tmp = MMG5_Get_basename(psl->namein);
1950 if ( metricData ) {
1951 fprintf(inm,"\"%s:metric\"\n",tmp);
1952 }
1953 else {
1954 fprintf(inm,"\"%s\"\n",tmp);
1955 }
1956 free(tmp); tmp = 0;
1957 }
1958 else {
1959 if ( metricData ) {
1960 fprintf(inm,"\"solution:metric\"\n");
1961 }
1962 else {
1963 fprintf(inm,"\"solution\"\n");
1964 }
1965 }
1966
1967 /* One real tag unused */
1968 fprintf(inm,"1\n");
1969 fprintf(inm,"0.0\n"); // time value: unused
1970
1971 /* Three integer tags */
1972 fprintf(inm,"3\n");
1973 fprintf(inm,"0\n"); // Time step: unused
1974 fprintf(inm,"%d\n",typ);
1975 fprintf(inm,"%" MMG5_PRId "\n",np);
1976
1979 if ( psl->size!= (psl->dim*(psl->dim+1))/2 ) {
1980 for (k=1; k<=mesh->np; k++) {
1981 ppt = &mesh->point[k];
1982 if ( !MG_VOK(ppt) ) continue;
1983
1984 iadr = k*psl->size;
1985 if ( psl->dim == 2 )
1986 dbuf[2] = 0; // z-component for a vector field
1987
1988 for ( i=0; i<psl->size; ++i )
1989 dbuf[i] = psl->m[iadr+i];
1990
1991 if ( !bin ) {
1992 fprintf(inm,"%" MMG5_PRId "",ppt->tmp);
1993 for ( i=0; i<typ; ++i )
1994 fprintf(inm," %lg",dbuf[i]);
1995 fprintf(inm,"\n");
1996 }
1997 else {
1998 fwrite(&ppt->tmp,MMG5_SW,1,inm);
1999 fwrite(&dbuf,MMG5_SD,typ,inm);
2000 }
2001 }
2002 }
2003 else {
2004 for (k=1; k<=mesh->np; k++) {
2005 ppt = &mesh->point[k];
2006 if ( !MG_VOK(ppt) ) continue;
2007
2008 if ( psl->dim == 3 ) {
2009 if ( metricData ) {
2010 assert(!mesh->nsols);
2011 MMG5_build3DMetric(mesh,psl,k,dbuf);
2012 }
2013 else {
2014 for (i=0; i<psl->size; i++) dbuf[i] = psl->m[psl->size*k+i];
2015 }
2016 }
2017
2018 if(!bin) {
2019 fprintf(inm,"%" MMG5_PRId "",ppt->tmp);
2020 if ( psl->dim==2 ) {
2021 iadr = k*psl->size;
2022 fprintf(inm," %.15lg %.15lg %.15lg %.15lg %.15lg %.15lg"
2023 " %.15lg %.15lg %.15lg \n",
2024 psl->m[iadr],psl->m[iadr+1],0.,psl->m[iadr+1],psl->m[iadr+2],0.,0.,0.,1.);
2025 }
2026 else {
2027 fprintf(inm," %.15lg %.15lg %.15lg %.15lg %.15lg %.15lg"
2028 " %.15lg %.15lg %.15lg \n", dbuf[0],dbuf[1],dbuf[2],
2029 dbuf[1],dbuf[3],dbuf[4],dbuf[2],dbuf[4],dbuf[5]);
2030 }
2031 }
2032 else {
2033 fwrite(&ppt->tmp,MMG5_SW,1,inm);
2034 if ( psl->dim==2 ) {
2035 iadr = k*psl->size;
2036 fwrite(&psl->m[iadr],MMG5_SD,2,inm);
2037 dbuf[0] = dbuf[1] = dbuf[2] = 0.;
2038 dbuf[3] = 1.;
2039 fwrite(&dbuf,MMG5_SD,1,inm);
2040 fwrite(&psl->m[iadr+1],MMG5_SD,2,inm);
2041 fwrite(&dbuf,MMG5_SD,4,inm);
2042 }
2043 else {
2044 fwrite(&dbuf[0],MMG5_SD,3,inm);
2045 fwrite(&dbuf[1],MMG5_SD,1,inm);
2046 fwrite(&dbuf[3],MMG5_SD,2,inm);
2047 fwrite(&dbuf[2],MMG5_SD,1,inm);
2048 fwrite(&dbuf[4],MMG5_SD,2,inm);
2049 }
2050 }
2051 }
2052 }
2053 if ( bin ) fprintf(inm,"\n");
2054 fprintf(inm,"$EndNodeData\n");
2055 }
2056 fclose(inm);
2057
2058 return 1;
2059}
2060
2080int MMG5_loadSolHeader( const char *filename,int meshDim,FILE **inm,int *ver,
2081 int *bin,int *iswp,MMG5_int *np,int *dim,int *nsols,int **type,
2082 long *posnp, int imprim) {
2083 int binch,bdim;
2084 int bpos,i;
2085 char *ptr,*data,chaine[MMG5_FILESTR_LGTH];
2086
2087 *posnp = 0;
2088 *bin = 0;
2089 *iswp = 0;
2090 *ver = 0;
2091
2092 MMG5_SAFE_CALLOC(data,strlen(filename)+6,char,return -1);
2093 strcpy(data,filename);
2094
2095 ptr = strstr(data,".mesh");
2096 if ( ptr ) *ptr = '\0';
2097
2098 ptr = strstr(data,".sol");
2099
2100 if ( !ptr ) {
2101 /* data contains the filename without extension */
2102 strcat(data,".solb");
2103 if (!(*inm = fopen(data,"rb")) ) {
2104 /* our file is not a .solb file, try with .sol ext */
2105 ptr = strstr(data,".solb");
2106 *ptr = '\0';
2107 strcat(data,".sol");
2108 if (!(*inm = fopen(data,"rb")) ) {
2109 if ( imprim >= 0 )
2110 fprintf(stderr," ** %s NOT FOUND. USE DEFAULT METRIC.\n",data);
2111 MMG5_SAFE_FREE(data);
2112 return 0;
2113 }
2114 } else {
2115 *bin = 1;
2116 }
2117 }
2118 else {
2119 ptr = strstr(data,".solb");
2120 if ( ptr ) *bin = 1;
2121
2122 if (!(*inm = fopen(data,"rb")) ) {
2123 if ( imprim >= 0 )
2124 fprintf(stderr," ** %s NOT FOUND. USE DEFAULT METRIC.\n",data);
2125 MMG5_SAFE_FREE(data);
2126 return 0;
2127 }
2128 }
2129 if ( imprim >= 0 )
2130 fprintf(stdout," %%%% %s OPENED\n",data);
2131 MMG5_SAFE_FREE(data);
2132
2133 /* read solution or metric */
2134 if(!*bin) {
2135 strcpy(chaine,"DDD");
2136 while(fscanf(*inm,"%127s",&chaine[0])!=EOF && strncmp(chaine,"End",strlen("End")) ) {
2137 if(!strncmp(chaine,"Dimension",strlen("Dimension"))) {
2138 MMG_FSCANF(*inm,"%d",dim);
2139 if ( *dim!=meshDim ) {
2140 fprintf(stderr,"BAD SOL DIMENSION: %d\n",*dim);
2141 fclose(*inm);
2142 return -1;
2143 }
2144 continue;
2145 } else if(!strncmp(chaine,"SolAtVertices",strlen("SolAtVertices"))) {
2146 MMG_FSCANF(*inm,"%" MMG5_PRId "",np);
2147 MMG_FSCANF(*inm,"%d",nsols);
2148 MMG5_SAFE_CALLOC(*type,*nsols,int,return -1);
2149 for ( i=0; i<*nsols; ++i ) {
2150 MMG_FSCANF(*inm,"%d",&((*type)[i]));
2151 }
2152 *posnp = ftell(*inm);
2153 break;
2154 }
2155 }
2156 } else {
2157 MMG_FREAD(&binch,MMG5_SW,1,*inm);
2158 if(binch==16777216) (*iswp)=1;
2159 else if(binch!=1) {
2160 fprintf(stdout,"BAD FILE ENCODING\n");
2161 fclose(*inm);
2162 return -1;
2163 }
2164 MMG_FREAD(ver,MMG5_SW,1,*inm);
2165 if ( *iswp ) *ver = MMG5_swapbin(*ver);
2166 while(fread(&binch,MMG5_SW,1,*inm)!=EOF && binch!=54 ) {
2167 if ( *iswp ) binch=MMG5_swapbin(binch);
2168 if(binch==54) break;
2169 if(binch==3) { //Dimension
2170 MMG_FREAD(&bdim,MMG5_SW,1,*inm); //NulPos=>20
2171 if ( *iswp ) bdim=MMG5_swapbin(bdim);
2172 MMG_FREAD(dim,MMG5_SW,1,*inm);
2173 if ( *iswp ) *dim=MMG5_swapbin(*dim);
2174 if ( *dim!=meshDim ) {
2175 fprintf(stderr,"BAD SOL DIMENSION: %d\n",*dim);
2176 printf(" Exit program.\n");
2177 fclose(*inm);
2178 return -1;
2179 }
2180 continue;
2181 } else if(binch==62) { //SolAtVertices
2182 MMG_FREAD(&binch,MMG5_SW,1,*inm); //Pos
2183 if ( *iswp ) binch=MMG5_swapbin(binch);
2184 MMG_FREAD(np,MMG5_SW,1,*inm);
2185 if ( *iswp ) *np=MMG5_swapbin(*np);
2186 MMG_FREAD(nsols,MMG5_SW,1,*inm); //nb sol
2187 if ( *iswp ) *nsols =MMG5_swapbin(*nsols);
2188
2189 MMG5_SAFE_CALLOC(*type,*nsols,int,return -1); //typSol
2190 for ( i=0; i<*nsols; ++i ) {
2191 MMG_FREAD(&((*type)[i]),MMG5_SW,1,*inm);
2192 if ( *iswp ) (*type)[i]=MMG5_swapbin((*type)[i]);
2193 }
2194 *posnp = ftell(*inm);
2195 break;
2196 } else {
2197 MMG_FREAD(&bpos,MMG5_SW,1,*inm); //Pos
2198 if ( *iswp ) bpos=MMG5_swapbin(bpos);
2199 rewind(*inm);
2200 fseek(*inm,bpos,SEEK_SET);
2201 }
2202 }
2203 }
2204
2205 return 1;
2206}
2207
2220int MMG5_readFloatSol3D(MMG5_pSol sol,FILE *inm,int bin,int iswp,int pos) {
2221 float fbuf[6],tmpf;
2222 int i;
2223
2224 switch ( sol->size ) {
2225 case 1: case 3:
2226 /* scalar or vector solution */
2227 for (i=0; i<sol->size; i++) {
2228 if(!bin){
2229 MMG_FSCANF(inm,"%f",&fbuf[0]);
2230 } else {
2231 MMG_FREAD(&fbuf[0],MMG5_SW,1,inm);
2232 if(iswp) fbuf[0]=MMG5_swapf(fbuf[0]);
2233 }
2234 sol->m[sol->size*pos+i] = fbuf[0];
2235 }
2236 break;
2237 case 6 :
2238 /* Tensor solution */
2239 if(!bin){
2240 for(i=0 ; i<sol->size ; i++)
2241 MMG_FSCANF(inm,"%f",&fbuf[i]);
2242 } else {
2243 for(i=0 ; i<sol->size ; i++) {
2244 MMG_FREAD(&fbuf[i],MMG5_SW,1,inm);
2245 if(iswp) fbuf[i]=MMG5_swapf(fbuf[i]);
2246 }
2247 }
2248 tmpf = fbuf[2];
2249 fbuf[2] = fbuf[3];
2250 fbuf[3] = tmpf;
2251 for (i=0; i<6; i++) sol->m[6*pos+i] = fbuf[i];
2252 break;
2253 }
2254 return 1;
2255}
2256
2269int MMG5_readDoubleSol3D(MMG5_pSol sol,FILE *inm,int bin,int iswp,MMG5_int pos) {
2270 double dbuf[6],tmpd;
2271 int i;
2272
2273 switch ( sol->size ) {
2274 case 1: case 3:
2275 /* scalar or vector solution */
2276 for (i=0; i<sol->size; i++) {
2277 if(!bin){
2278 MMG_FSCANF(inm,"%lf",&dbuf[i]);
2279 } else {
2280 MMG_FREAD(&dbuf[i],MMG5_SD,1,inm);
2281 if(iswp) dbuf[i]=MMG5_swapd(dbuf[i]);
2282 }
2283 sol->m[sol->size*pos+i] = dbuf[i];
2284 }
2285 break;
2286
2287 case 6 :
2288 /* tensor solution */
2289 if(!bin){
2290 for(i=0 ; i<sol->size ; i++)
2291 MMG_FSCANF(inm,"%lf",&dbuf[i]);
2292 } else {
2293 for(i=0 ; i<sol->size ; i++) {
2294 MMG_FREAD(&dbuf[i],MMG5_SD,1,inm);
2295 if(iswp) dbuf[i]=MMG5_swapd(dbuf[i]);
2296 }
2297 }
2298 tmpd = dbuf[2];
2299 dbuf[2] = dbuf[3];
2300 dbuf[3] = tmpd;
2301 for (i=0; i<sol->size; i++) sol->m[6*pos+i] = dbuf[i];
2302 break;
2303 }
2304 return 1;
2305}
2306
2319 MMG5_int pos,int metricData) {
2320 double dbuf[6],tmp;
2321 int i;
2322
2323 switch ( sol->size ) {
2324 case 1: case 3:
2325 /* scalar or vector solution */
2326 for (i=0; i<sol->size; i++) {
2327 for (i=0; i<sol->size; i++) dbuf[i] = sol->m[sol->size*pos+i];
2328 if(!bin){
2329 for (i=0; i<sol->size; i++)
2330 fprintf(inm," %.15lg",dbuf[i]);
2331 } else {
2332 for(i=0; i<sol->size; i++)
2333 fwrite((unsigned char*)&dbuf[i],MMG5_SD,1,inm);
2334 }
2335 }
2336 break;
2337
2338 case 6 :
2339 /* tensor solution */
2340 if ( metricData )
2341 {
2342 MMG5_build3DMetric(mesh,sol,pos,dbuf);
2343 }
2344 else
2345 {
2346 for (i=0; i<sol->size; i++)
2347 {
2348 dbuf[i] = sol->m[6*pos+i];
2349 }
2350 }
2351
2352 tmp = dbuf[2];
2353 dbuf[2] = dbuf[3];
2354 dbuf[3] = tmp;
2355
2356 if(!bin) {
2357 for(i=0; i<sol->size; i++)
2358 fprintf(inm," %.15lg",dbuf[i]);
2359 } else {
2360 for(i=0; i<sol->size; i++)
2361 fwrite((unsigned char*)&dbuf[i],MMG5_SD,1,inm);
2362 }
2363 break;
2364 }
2365}
2366
2388 FILE **inm,int ver,int *bin,MMG5_int *bpos,MMG5_int np,int dim,
2389 int nsols,int *entities,int *type,int *size) {
2390 MMG5_pPoint ppt;
2391 int binch;
2392 MMG5_int k;
2393 char *ptr,*data,chaine[MMG5_FILESTR_LGTH];
2394
2395 *bin = 0;
2396
2397 MMG5_SAFE_CALLOC(data,strlen(filename)+6,char,return 0);
2398 strcpy(data,filename);
2399 ptr = strstr(data,".sol");
2400 if ( ptr ) {
2401 // filename contains the solution extension
2402 ptr = strstr(data,".solb");
2403
2404 if ( ptr ) *bin = 1;
2405
2406 if( !(*inm = fopen(data,"wb")) ) {
2407 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
2408 MMG5_SAFE_FREE(data);
2409 return 0;
2410 }
2411 }
2412 else
2413 {
2414 // filename don't contains the solution extension
2415 ptr = strstr(data,".mesh");
2416 if ( ptr ) *ptr = '\0';
2417
2418 strcat(data,".sol");
2419 if (!(*inm = fopen(data,"wb")) ) {
2420 ptr = strstr(data,".solb");
2421 *ptr = '\0';
2422 strcat(data,".sol");
2423 if (!(*inm = fopen(data,"wb")) ) {
2424 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
2425 MMG5_SAFE_FREE(data);
2426 return 0;
2427 }
2428 else *bin = 1;
2429 }
2430 }
2431
2432 if ( mesh->info.imprim >= 0 )
2433 fprintf(stdout," %%%% %s OPENED\n",data);
2434 MMG5_SAFE_FREE(data);
2435
2436 /*entete fichier*/
2437 binch=(*bpos)=0;
2438 if(!*bin) {
2439 strcpy(&chaine[0],"MeshVersionFormatted\n");
2440 fprintf(*inm,"%s %d",chaine,ver);
2441 strcpy(&chaine[0],"\n\nDimension\n");
2442 fprintf(*inm,"%s %d",chaine,dim);
2443 } else {
2444 binch = 1; //MeshVersionFormatted
2445 fwrite(&binch,MMG5_SW,1,*inm);
2446 binch = ver; //version
2447 fwrite(&binch,MMG5_SW,1,*inm);
2448 binch = 3; //Dimension
2449 fwrite(&binch,MMG5_SW,1,*inm);
2450 (*bpos) = 20; //Pos
2451 fwrite(bpos,MMG5_SW,1,*inm);
2452 binch = dim;
2453 fwrite(&binch,MMG5_SW,1,*inm);
2454 }
2455
2456 np = 0;
2457 for (k=1; k<=mesh->np; k++) {
2458 ppt = &mesh->point[k];
2459 if ( MG_VOK(ppt) ) np++;
2460 }
2461
2462 /* Count the number of solutions at vertices */
2463 int npointSols = 0;
2464 for (k=0; k<nsols; ++k ) {
2465 if ( (entities[k]==MMG5_Noentity) || (entities[k]==MMG5_Vertex) ) {
2466 ++npointSols;
2467 }
2468 }
2469
2470 /* Sol at vertices header */
2471 if(!*bin) {
2472 strcpy(&chaine[0],"\n\nSolAtVertices\n");
2473 fprintf(*inm,"%s",chaine);
2474 fprintf(*inm,"%" MMG5_PRId "\n",np);
2475 fprintf(*inm,"%d",npointSols);
2476 for (k=0; k<nsols; ++k ) {
2477 if ( (entities[k]!=MMG5_Noentity) && (entities[k]!=MMG5_Vertex) ) {
2478 continue;
2479 }
2480 fprintf(*inm," %d",type[k]);
2481 }
2482 fprintf(*inm,"\n");
2483 } else {
2484 binch = 62; //Vertices
2485 fwrite(&binch,MMG5_SW,1,*inm);
2486 (*bpos) += 16;
2487
2488 for (k=0; k<nsols; ++k ) {
2489 if ( (entities[k]!=MMG5_Noentity) && (entities[k]!=MMG5_Vertex) ) {
2490 continue;
2491 }
2492 (*bpos) += 4 + (size[k]*ver)*4*np; //Pos
2493 }
2494 fwrite(bpos,MMG5_SW,1,*inm);
2495
2496 fwrite(&np,MMG5_SW,1,*inm);
2497 fwrite(&npointSols,MMG5_SW,1,*inm);
2498 for (k=0; k<nsols; ++k ) {
2499 if ( (entities[k]!=MMG5_Noentity) && (entities[k]!=MMG5_Vertex) ) {
2500 continue;
2501 }
2502 fwrite(&type[k],MMG5_SW,1,*inm);
2503 }
2504 }
2505
2506 return 1;
2507}
2508
2525 FILE *inm,int ver,int bin,MMG5_int *bpos,
2526 int nsols,int nsolsAtTriangles,
2527 int *entities,int *type,int *size) {
2528 MMG5_pTria pt;
2529 int binch;
2530 MMG5_int k,nt;
2531
2532 nt = 0;
2533 for (k=1; k<=mesh->nt; k++) {
2534 pt = &mesh->tria[k];
2535 if ( MG_EOK(pt) ) {
2536 ++nt;
2537 }
2538 }
2539
2540 /* Sol at vertices header */
2541 if(!bin) {
2542 fprintf(inm,"\n\nSolAtTriangles\n");
2543 fprintf(inm,"%" MMG5_PRId "\n",nt);
2544 fprintf(inm,"%d",nsolsAtTriangles);
2545 for (k=0; k<nsols; ++k ) {
2546 if ( entities[k]!=MMG5_Triangle ) {
2547 continue;
2548 }
2549 fprintf(inm," %d",type[k]);
2550 }
2551 fprintf(inm,"\n");
2552 } else {
2553 binch = 64; //SolsAtTriangles
2554 fwrite(&binch,MMG5_SW,1,inm);
2555 (*bpos) += 16;
2556
2557 for (k=0; k<nsols; ++k ) {
2558 if ( entities[k]!=MMG5_Triangle ) {
2559 continue;
2560 }
2561 (*bpos) += 4 + (size[k]*ver)*4*nt; //Pos
2562 }
2563 fwrite(bpos,MMG5_SW,1,inm);
2564
2565 fwrite(&nt,MMG5_SW,1,inm);
2566 fwrite(&nsolsAtTriangles,MMG5_SW,1,inm);
2567 for (k=0; k<nsols; ++k ) {
2568 if ( entities[k]!=MMG5_Triangle ) {
2569 continue;
2570 }
2571 fwrite(&type[k],MMG5_SW,1,inm);
2572 }
2573 }
2574
2575 return 1;
2576}
2577
2594 FILE *inm,int ver,int bin,MMG5_int *bpos,
2595 int nsols,int nsolsAtTetra,
2596 int *entities,int *type,int *size) {
2597 MMG5_pTetra pt;
2598 int binch;
2599 MMG5_int k,ne;
2600
2601 ne = 0;
2602 for (k=1; k<=mesh->ne; k++) {
2603 pt = &mesh->tetra[k];
2604 if ( MG_EOK(pt) ) {
2605 ++ne;
2606 }
2607 }
2608
2609 /* Sol at vertices header */
2610 if(!bin) {
2611 fprintf(inm,"\n\nSolAtTetrahedra\n");
2612 fprintf(inm,"%" MMG5_PRId "\n",ne);
2613 fprintf(inm,"%d",nsolsAtTetra);
2614 for (k=0; k<nsols; ++k ) {
2615 if ( entities[k]!=MMG5_Tetrahedron ) {
2616 continue;
2617 }
2618 fprintf(inm," %d",type[k]);
2619 }
2620 fprintf(inm,"\n");
2621 } else {
2622 binch = 66; //SolsAtTetrahedron
2623 fwrite(&binch,MMG5_SW,1,inm);
2624 (*bpos) += 16;
2625
2626 for (k=0; k<nsols; ++k ) {
2627 if ( entities[k]!=MMG5_Tetrahedron ) {
2628 continue;
2629 }
2630 (*bpos) += 4 + (size[k]*ver)*4*ne; //Pos
2631 }
2632 fwrite(bpos,MMG5_SW,1,inm);
2633
2634 fwrite(&ne,MMG5_SW,1,inm);
2635 fwrite(&nsolsAtTetra,MMG5_SW,1,inm);
2636 for (k=0; k<nsols; ++k ) {
2637 if ( entities[k]!=MMG5_Tetrahedron ) {
2638 continue;
2639 }
2640 fwrite(&type[k],MMG5_SW,1,inm);
2641 }
2642 }
2643
2644 return 1;
2645}
2646
2661int MMG5_chkMetricType(MMG5_pMesh mesh,int *type, int *entities, FILE *inm) {
2662
2663 /* Metric can be provided only on vertices */
2664 if ( (*entities != MMG5_Vertex) && (*entities != MMG5_Noentity) ) {
2665 fprintf(stderr," ## Error: %s: Metric should apply on vertices.\n"
2666 " If your input file is at a non Medit"
2667 " file format, please ensure to remove non metric fields from your"
2668 " file and that the metric field"
2669 " contains the \":metric\" string.\n",__FILE__);
2670 if ( inm ) fclose(inm);
2671 return -1;
2672 }
2673
2674 /* 1: scalar solution (isotropic metric or ls function,
2675 2: vector field (displacement in Lagrangian mode),
2676 3: anisotropic metric */
2677 if ( mesh->info.lag == -1 ) {
2678 if ( type[0]!=1 && type[0]!=3) {
2679 fprintf(stderr," ** DATA TYPE IGNORED %d \n",type[0]);
2680 fprintf(stderr," ## Error: %s: if your input file is at a non Medit"
2681 " file format, please ensure that the metric field"
2682 " contains the \":metric\" string.\n",__FILE__);
2683 if ( inm ) fclose(inm);
2684 return -1;
2685 }
2686 }
2687 else {
2688 if ( type[0] != 2 ) {
2689 fprintf(stderr," ** MISMATCH DATA TYPE FOR LAGRANGIAN MODE %d \n",
2690 type[0]);
2691 if ( inm ) fclose(inm);
2692 return -1;
2693 }
2694 }
2695 return 1;
2696}
2697
2706 if ( abs(mesh->info.imprim) > 3 ) {
2707 if ( met->size == 1 )
2708 fprintf(stdout," NUMBER OF SCALAR VALUES %8" MMG5_PRId "\n",met->np);
2709 else if ( met->size == 3 )
2710 fprintf(stdout," NUMBER OF VECTOR VALUES %8" MMG5_PRId "\n",met->np);
2711 else
2712 fprintf(stdout," NUMBER OF TENSOR VALUES %8" MMG5_PRId "\n",met->np);
2713 }
2714}
2715
2724 int j;
2725
2726 if ( abs(mesh->info.imprim) > 3 ) {
2727 fprintf(stdout," NUMBER OF SOLUTIONS PER ENTITY %8d\n",mesh->nsols);
2728 fprintf(stdout," TYPE OF SOLUTIONS:\n ");
2729 for ( j=0; j<mesh->nsols; ++j ) {
2730 if ( (*sol)[j].size == 1 )
2731 fprintf(stdout," SCALAR");
2732 else if ( (*sol)[j].size == 3 )
2733 fprintf(stdout," VECTOR");
2734 else
2735 fprintf(stdout," TENSOR");
2736 }
2737 fprintf(stdout,"\n");
2738 }
2739}
2740
2742 FILE* inm;
2743 MMG5_pPoint ppt;
2744 int i;
2745 MMG5_int k,np;
2746 char *ptr,*data;
2747
2748 if ( !mesh->np ) {
2749 return 1;
2750 }
2751
2752 if ( (!filename) || !(*filename) ) {
2754 }
2755 if ( (!filename) || !(*filename) ) {
2756 printf("\n ## Error: %s: unable to save a file without a valid filename\n.",
2757 __func__);
2758 return 0;
2759 }
2760
2761 /* Name of file */
2762 MMG5_SAFE_CALLOC(data,strlen(filename)+7,char,return 0);
2763 strcpy(data,filename);
2764 ptr = strstr(data,".node");
2765 if ( ptr ) {
2766 *ptr = '\0';
2767 }
2768
2769 /* Add .node ext */
2770 strcat(data,".node");
2771 if( !(inm = fopen(data,"wb")) ) {
2772 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
2773 MMG5_SAFE_FREE(data);
2774 return 0;
2775 }
2776
2777 fprintf(stdout," %%%% %s OPENED\n",data);
2778 MMG5_SAFE_FREE(data);
2779
2780 /* Write vertices */
2781 np = 0;
2782 for (k=1; k<=mesh->np; k++) {
2783 ppt = &mesh->point[k];
2784 if ( MG_VOK(ppt) ) {
2785 ++np;
2786 ppt->tmp = np;
2787 }
2788 }
2789
2790 /* Save node number, dim, no attributes, 1 bdy marker */
2791 fprintf(inm, "%" MMG5_PRId " %d %d %d\n\n",np,mesh->dim,0,1);
2792
2793 for ( k=1; k<=mesh->np; ++k ) {
2794 /* Save node idx */
2795 ppt = &mesh->point[k];
2796 if ( MG_VOK(ppt) ) {
2797 fprintf(inm, "%" MMG5_PRId " ",ppt->tmp);
2798
2799 /* Save coordinates */
2800 for ( i=0; i<mesh->dim; ++i ) {
2801 fprintf(inm, " %.15lf",ppt->c[i]);
2802 }
2803
2804 /* Save bdy marker */
2805 fprintf(inm, " %" MMG5_PRId "\n",ppt->ref);
2806 }
2807 }
2808
2809 fprintf(stdout," NUMBER OF VERTICES %8" MMG5_PRId "\n",np);
2810
2811 fclose(inm);
2812
2813 return 1;
2814}
2815
2816int MMG5_saveEdge(MMG5_pMesh mesh,const char *filename,const char *ext) {
2817 FILE* inm;
2818 MMG5_pEdge pt;
2819 size_t na_tot;
2820 int polyfile;
2821 MMG5_int k;
2822 char *ptr_c = (char*)mesh->edge;
2823 char *ptr,*data;
2824
2825 if ( !mesh->edge ) {
2826 return 1;
2827 }
2828 if ( !mesh->na ) {
2829 return 1;
2830 }
2831
2832 if ( (!filename) || !(*filename) ) {
2834 }
2835 if ( (!filename) || !(*filename) ) {
2836 printf("\n ## Error: %s: unable to save a file without a valid filename\n.",
2837 __func__);
2838 return 0;
2839 }
2840
2841 /* Name of file */
2842 MMG5_SAFE_CALLOC(data,strlen(filename)+strlen(ext),char,return 0);
2843 strcpy(data,filename);
2844 ptr = strstr(data,".node");
2845 if ( ptr ) {
2846 *ptr = '\0';
2847 }
2848
2849 /* Add file ext */
2850 strcat(data,ext);
2851 if( !(inm = fopen(data,"wb")) ) {
2852 fprintf(stderr," ** UNABLE TO OPEN %s.\n",data);
2853 MMG5_SAFE_FREE(data);
2854 return 0;
2855 }
2856
2857 fprintf(stdout," %%%% %s OPENED\n",data);
2858 MMG5_SAFE_FREE(data);
2859
2860 /* For .poly file, add header */
2861 if ( !strcmp(ext,".poly") ) {
2862 polyfile = 1;
2863 }
2864 else {
2865 polyfile = 0;
2866 }
2867
2868 if ( polyfile ) {
2869 /* Save 0 nodes (saved in a separated .node file), dim, 0 attributes, 1 bdy
2870 * marker */
2871 fprintf(inm, "0 %d 0 1\n",mesh->dim);
2872 }
2873
2874 /* Get either the number of boundary edges or the total number of edges
2875 * (depending if they have been append to the bdy edges, if yes, edges 1->na
2876 * are bdy, na->na_tot are internal. */
2877
2878 /* Get size of the array in octets */
2879 ptr_c = ptr_c-sizeof(size_t);
2880 na_tot = (*((size_t*)ptr_c));
2881 /* Recover number of edges allocated */
2882 na_tot /= sizeof(MMG5_Edge);
2883 /* Array is allocated at size na+1, recover na */
2884 --na_tot;
2885
2886 /* Save node number, dim, no attributes, 1 bdy marker */
2887 fprintf(inm, "%zu %d\n",na_tot,1);
2888
2889 for ( k=1; k<=na_tot; ++k ) {
2890 /* Save edge idx */
2891 fprintf(inm, "%" MMG5_PRId " ",k);
2892
2893 pt = &mesh->edge[k];
2894
2895 /* Save connectivity */
2896 fprintf(inm,"%" MMG5_PRId " %" MMG5_PRId " %" MMG5_PRId "\n",mesh->point[pt->a].tmp,mesh->point[pt->b].tmp,pt->ref);
2897
2898 }
2899
2900 /* For .poly file, add last line: 0 holes */
2901 if ( polyfile ) {
2902 fprintf(inm, "0 \n");
2903 }
2904
2905 fprintf(stdout," NUMBER OF EDGES %8zu\n",na_tot);
2906
2907 fclose(inm);
2908
2909 return 1;
2910}
char * MMG5_Get_basename(char *path)
int MMG5_Set_inputSolName(MMG5_pMesh mesh, MMG5_pSol sol, const char *solin)
int ier
tmp[*strlen0]
MMG5_pMesh MMG5_pSol * sol
if(!ier) exit(EXIT_FAILURE)
MMG5_pMesh * mesh
MMG5_pMesh char * filename
int 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_build3DMetric(MMG5_pMesh mesh, MMG5_pSol sol, MMG5_int ip, double dbuf[6])
Definition: inout.c:1528
void MMG5_printSolStats(MMG5_pMesh mesh, MMG5_pSol *sol)
Definition: inout.c:2723
float MMG5_swapf(float sbin)
Definition: inout.c:80
int MMG5_check_readedMesh(MMG5_pMesh mesh, MMG5_int nref)
Definition: inout.c:526
static int MMG5_countBinaryElts(FILE **inm, const int nelts, const int iswp, int *np, int *na, int *nt, int *nq, int *ne, int *npr)
Definition: inout.c:127
int MMG5_saveSolAtTrianglesHeader(MMG5_pMesh mesh, FILE *inm, int ver, int bin, MMG5_int *bpos, int nsols, int nsolsAtTriangles, int *entities, int *type, int *size)
Definition: inout.c:2524
int MMG5_saveNode(MMG5_pMesh mesh, const char *filename)
Definition: inout.c:2741
int MMG5_saveEdge(MMG5_pMesh mesh, const char *filename, const char *ext)
Definition: inout.c:2816
int MMG5_saveSolHeader(MMG5_pMesh mesh, const char *filename, FILE **inm, int ver, int *bin, MMG5_int *bpos, MMG5_int np, int dim, int nsols, int *entities, int *type, int *size)
Definition: inout.c:2387
int MMG5_loadMshMesh_part1(MMG5_pMesh mesh, const char *filename, FILE **inm, long *posNodes, long *posElts, long **posNodeData, int *bin, int *iswp, MMG5_int *nelts, int *nsols)
Definition: inout.c:278
void MMG5_writeDoubleSol3D(MMG5_pMesh mesh, MMG5_pSol sol, FILE *inm, int bin, MMG5_int pos, int metricData)
Definition: inout.c:2318
MMG5_int MMG5_swapbin_int(MMG5_int sbin)
Definition: inout.c:61
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:2593
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
@ MMG5_Noentity
Definition: libmmgtypes.h:223
@ MMG5_Vertex
Definition: libmmgtypes.h:224
@ MMG5_Tetrahedron
Definition: libmmgtypes.h:227
@ MMG5_Triangle
Definition: libmmgtypes.h:226
#define MG_GEO
#define MMG5_SAFE_CALLOC(ptr, size, type, law)
double MMG2D_quickarea(double a[2], double b[2], double c[2])
Definition: tools.c:971
#define MG_EOK(pt)
#define MG_NUL
#define MMG5_SD
#define MMG5_ADD_MEM(mesh, size, message, law)
#define MG_SIN(tag)
#define MMG_FREAD(ptr, size, count, stream)
#define MMG5_SW
#define MG_VOK(ppt)
double MMG5_orvol(MMG5_pPoint point, MMG5_int *v)
Definition: tools.c:951
#define MG_NOSURF
#define MG_NOM
#define MMG5_FILESTR_LGTH
#define MG_REF
#define MMG5_SAFE_FREE(ptr)
#define MMG_FSCANF(stream, format,...)
#define MMG5_DEL_MEM(mesh, ptr)
#define MMG5_SAFE_RECALLOC(ptr, prevSize, newSize, type, message, law)
Structure to store edges of a MMG mesh.
Definition: libmmgtypes.h:305
MMG5_int b
Definition: libmmgtypes.h:306
MMG5_int ref
Definition: libmmgtypes.h:307
int16_t tag
Definition: libmmgtypes.h:310
MMG5_int a
Definition: libmmgtypes.h:306
int8_t iso
Definition: libmmgtypes.h:534
MMG5_int isoref
Definition: libmmgtypes.h:521
uint8_t metRidTyp
Definition: libmmgtypes.h:547
int8_t lag
Definition: libmmgtypes.h:540
MMG mesh structure.
Definition: libmmgtypes.h:605
MMG5_pQuad quadra
Definition: libmmgtypes.h:648
MMG5_Info info
Definition: libmmgtypes.h:651
MMG5_int xt
Definition: libmmgtypes.h:620
MMG5_int ne
Definition: libmmgtypes.h:612
MMG5_pPoint point
Definition: libmmgtypes.h:641
char * nameout
Definition: libmmgtypes.h:653
MMG5_pPrism prism
Definition: libmmgtypes.h:645
MMG5_int nquad
Definition: libmmgtypes.h:613
MMG5_int npmax
Definition: libmmgtypes.h:612
MMG5_pxPoint xpoint
Definition: libmmgtypes.h:642
MMG5_int xp
Definition: libmmgtypes.h:620
MMG5_int nei
Definition: libmmgtypes.h:612
MMG5_pTetra tetra
Definition: libmmgtypes.h:643
MMG5_int nt
Definition: libmmgtypes.h:612
MMG5_pTria tria
Definition: libmmgtypes.h:647
MMG5_int np
Definition: libmmgtypes.h:612
MMG5_pEdge edge
Definition: libmmgtypes.h:649
MMG5_int nti
Definition: libmmgtypes.h:612
MMG5_int npi
Definition: libmmgtypes.h:612
MMG5_int nai
Definition: libmmgtypes.h:612
MMG5_int nprism
Definition: libmmgtypes.h:613
MMG5_int na
Definition: libmmgtypes.h:612
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
MMG5_int v[6]
Definition: libmmgtypes.h:463
MMG5_int ref
Definition: libmmgtypes.h:464
MMG5_int ref
Definition: libmmgtypes.h:368
MMG5_int v[4]
Definition: libmmgtypes.h:367
MMG5_int npi
Definition: libmmgtypes.h:667
int entities
Definition: libmmgtypes.h:670
MMG5_int npmax
Definition: libmmgtypes.h:666
char * namein
Definition: libmmgtypes.h:673
double * m
Definition: libmmgtypes.h:671
MMG5_int np
Definition: libmmgtypes.h:665
MMG5_int v[4]
Definition: libmmgtypes.h:403
MMG5_int ref
Definition: libmmgtypes.h:404
MMG5_int edg[3]
Definition: libmmgtypes.h:339
MMG5_int ref
Definition: libmmgtypes.h:335
MMG5_int v[3]
Definition: libmmgtypes.h:334
double n1[3]
Definition: libmmgtypes.h:295