File: MAdOutput.cc

package info (click to toggle)
madlib 1.3.0-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 8,196 kB
  • sloc: cpp: 39,851; sh: 10,041; makefile: 473
file content (816 lines) | stat: -rw-r--r-- 26,970 bytes parent folder | download | duplicates (6)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
// -------------------------------------------------------------------
// MAdLib - Copyright (C) 2008-2009 Universite catholique de Louvain
//
// See the Copyright.txt and License.txt files for license information. 
// You should have received a copy of these files along with MAdLib. 
// If not, see <http://www.madlib.be/license/>
//
// Please report all bugs and problems to <contrib@madlib.be>
//
// Authors: Gaetan Compere, Jean-Francois Remacle
// -------------------------------------------------------------------

#include "MAdOutput.h"
#include "MeanRatioEvaluator.h"
#include "OrientedMeanRatioEvaluator.h"
#include "MAdMessage.h"
#include "MathUtils.h"
#include "LocalSizeField.h"
#include "DistanceFunction.h"
#include "MAdResourceManager.h"
#include "MeshSizeBase.h"

#include <iostream>
using std::cout;
using std::cerr;
using std::endl;
#include <vector>
using std::vector;
using std::string;
using std::set;

namespace MAd {

  // ----------------------------------------------------------------------
  enum dataType {
    SCALAR,
    VECTORIAL
  };

  dataType getOutputDataType(MAdOutputData t) {
    switch (t) {
    case OD_CONSTANT:
    case OD_MEANRATIO:
    case OD_ORIENTEDMEANRATIO:
    case OD_SIZEFIELD_MEAN:
    case OD_SIZEFIELD_MIN:
    case OD_SIZEFIELD_MAX:
    case OD_DIMENSION:
    case OD_ITERORDER:
    case OD_CURVATURE_DIV:
    case OD_CURVATURE_MAX:
    case OD_CURVATURE_MIN:
    case OD_DISTANCE:
      { return SCALAR; break; }
    case OD_CURVATURE_MAX_VEC:
    case OD_CURVATURE_MIN_VEC:
    case OD_ANISO_SF_AXIS0:
    case OD_ANISO_SF_AXIS1:
    case OD_ANISO_SF_AXIS2:
      { return VECTORIAL; break; }
    }
    throw;
  };

  std::string getOutputName( MAdOutputData type ) {
    switch (type) {
    case OD_CONSTANT: return "constant field";
    case OD_MEANRATIO: return "mean ratio";
    case OD_ORIENTEDMEANRATIO: return "oriented mean ratio";
    case OD_SIZEFIELD_MEAN: return "mean size (among the 3 directions) in the size field";
    case OD_SIZEFIELD_MIN: return "minimum size (among the 3 directions) in the size field";
    case OD_SIZEFIELD_MAX: return "maximum size (among the 3 directions) in the size field";
    case OD_DIMENSION: return "element dimension";
    case OD_ITERORDER: return "element id regarding place in iterator";
    case OD_CURVATURE_DIV: return "divergence of the curvature";
    case OD_CURVATURE_MAX: return "surface maximum curvature (scalar field)";
    case OD_CURVATURE_MIN: return "surface minimum curvature (scalar field)";
    case OD_CURVATURE_MAX_VEC: return "surface maximum curvature (vectorial field)";
    case OD_CURVATURE_MIN_VEC: return "surface minimum curvature (vectorial field)";
    case OD_ANISO_SF_AXIS0: return "anisotropic size field along first direction";
    case OD_ANISO_SF_AXIS1: return "anisotropic size field along second direction";
    case OD_ANISO_SF_AXIS2: return "anisotropic size field along third direction";
    case OD_DISTANCE: return "distance to the walls";
    }
    return "unknown output type";
  };

  // ----------------------------------------------------------------------
  double* getData(MAdOutputData type, const pEntity pe, const pSField sf=NULL, int id=0)
  {
    int dim = EN_type(pe);
    double * result;

    switch (type) {

    case OD_DIMENSION: {
      result = new double[4];
      for (int i=0; i<4; i++) result[i] = (double)dim;
      return result;
    }

    case OD_CONSTANT: {
      result = new double[4];
      for (int i=0; i<4; i++) result[i] = 1.;
      return result;
    }

    case OD_MEANRATIO: {
      if ( !sf ) {
        MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                    "No size field given to compute mean ratio");
      }
      if ( sf->getType() != DISCRETESFIELD ) {
        MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                    "Discrete size field required for quality evaluation");
      }
      meanRatioEvaluator evalu((DiscreteSF*)sf);
      double tmp;
      switch (dim) {
      case 3:
        result = new double[4];
        evalu.R_shape((pRegion)pe,&tmp);
        for (int i=0; i<4; i++) result[i] = pow(tmp,1./3.);
        return result;
      case 2:
        result = new double[3];
        evalu.F_shape((pFace)pe,0,&tmp); 
        for (int i=0; i<3; i++) result[i] = sqrt(tmp);
        return result;
      default:
        cerr << "Error in outputs: getData on an element if dimension inferior to 2\n";
        throw;
      }
    }

    case OD_ORIENTEDMEANRATIO: {
      if ( !sf ) {
        MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                    "No size field given to compute mean ratio");
      }
      if ( sf->getType() != DISCRETESFIELD ) {
        MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                    "Discrete size field required for quality evaluation");
      }
      orientedMeanRatioEvaluator evalu((DiscreteSF*)sf);
      double tmp;
      switch (dim) {
      case 3:
        result = new double[4];
        evalu.R_shape((pRegion)pe,&tmp);
        for (int i=0; i<4; i++) result[i] = tmp;
        return result;
      case 2:
        result = new double[3];
        evalu.F_shape((pFace)pe,0,&tmp); 
        for (int i=0; i<3; i++) result[i] = tmp;
        return result;
      default:
        cerr << "Error in outputs: getData on an element if dimension inferior to 2\n";
        throw;
      }
    }

    case OD_ITERORDER: {
      result = new double[4];
      for (int i=0; i<4; i++) result[i] = (double)id;
      return result;
    }

    case OD_SIZEFIELD_MEAN:
    case OD_SIZEFIELD_MIN:
    case OD_SIZEFIELD_MAX: {
      if ( !sf ) {
        MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                    "No size field given");
      }
      void * temp; 
      pPList verts;
      int k;
      pVertex pv;
      pMSize size;
      double h;
      switch (dim) {
      case 3:
        result = new double[4];
        verts = R_vertices((pRegion)pe);
        temp = 0 ; k=0;
        while ( ( pv = (pVertex)PList_next(verts,&temp) ) )
          {
            size = sf->getSize(pv);
            if ( type == OD_SIZEFIELD_MEAN ) h = size->getMeanLength();
            if ( type == OD_SIZEFIELD_MIN )  h = size->getMinLength();
            if ( type == OD_SIZEFIELD_MAX )  h = size->getMaxLength();
            if (size) delete size;
            result[k] = h;
            k++;
          }
        PList_delete(verts);
        return result;
      case 2:
        result = new double[3];
        verts = F_vertices( (pFace)pe, 1);
        temp = 0 ; k=0;
        while ( ( pv = (pVertex)PList_next(verts,&temp) ) )
          {
            size = sf->getSize(pv);
            if ( type == OD_SIZEFIELD_MEAN ) h = size->getMeanLength();
            if ( type == OD_SIZEFIELD_MIN )  h = size->getMinLength();
            if ( type == OD_SIZEFIELD_MAX )  h = size->getMaxLength();
            if (size) delete size;
            result[k] = h;
            k++;
          }
        PList_delete(verts);
        return result;
      default:
        cerr << "Error in outputs: getData on an element if dimension inferior to 2\n";
        throw;
      }
    }

    case OD_CURVATURE_DIV: {
      if ( dim < 3 )
        {
#ifdef _HAVE_GMSH_
          result = new double[4];
          pGEntity pge = EN_whatIn(pe);
          if ( GEN_type(pge) != 2 ) {
            for (int i=0; i<4; i++) result[i] = -1.;
            return result;
          }
          double u[4][2];
          F_params((pFace)pe,u);
          for (int iV=0; iV<F_numVertices((pFace)pe); iV++)
            {
              result[iV] = GF_curvatureDiv((pGFace)pge, u[iV], 1.e6);
            }
          return result;
#else
          MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                      "Gmsh required for divergence of surface curvature");
#endif
        }
      else
        {
          if ( !sf || ( sf->getType() != LOCALSFIELD ) ) {
            MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                        "A local size field is required to output divergence of the curvature on elements");
          }
          result = new double[4];
          LocalSizeField * sf_cast = (LocalSizeField *) sf;
          pRegion pr = (pRegion) pe;
          for (int iV=0; iV<4; iV++)
            {
              if ( sf_cast->getCurvature(R_vertex(pr,iV), &(result[iV])) ) {}
              else result[iV] = -1.;
            }
          return result;
        }
      return NULL;
    }

    case OD_CURVATURE_MAX: {
#ifdef _HAVE_GMSH_
      if ( dim != 2 ) {
        MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                    "Cannot output curvatures on other entities than faces");
      }
      result = new double[4];
      pGEntity pge = EN_whatIn(pe);
      if ( GEN_type(pge) != 2 ) {
        for (int i=0; i<4; i++) result[i] = -1.;
        return result;
      }
      double u[4][2];
      F_params((pFace)pe,u);
      for (int iV=0; iV<F_numVertices((pFace)pe); iV++)
        {
          double dirMax[3], dirMin[3], curvMax, curvMin;
          GF_curvatures((pGFace)pge, u[iV],
                        dirMax, dirMin, &curvMax, &curvMin, 1.e6);
          result[iV] = curvMax;
        }
      return result;
#else
      MAdMsgSgl::instance().error(__LINE__,__FILE__,"Gmsh required for curvature");
#endif
    }

    case OD_CURVATURE_MIN: {
#ifdef _HAVE_GMSH_
      if ( dim != 2 ) {
        MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                    "Cannot output curvatures on other entities than faces");
      }
      result = new double[4];
      pGEntity pge = EN_whatIn(pe);
      if ( GEN_type(pge) != 2 ) {
        for (int i=0; i<4; i++) result[i] = -1.;
        return result;
      }
      double u[4][2];
      F_params((pFace)pe,u);
      for (int iV=0; iV<F_numVertices((pFace)pe); iV++)
        {
          double dirMax[3], dirMin[3], curvMax, curvMin;
          GF_curvatures((pGFace)pge, u[iV],
                        dirMax, dirMin, &curvMax, &curvMin, 1.e6);
          result[iV] = curvMin;
        }
      return result;
#else
      MAdMsgSgl::instance().error(__LINE__,__FILE__,"Gmsh required for curvature");
#endif
    }

    case OD_CURVATURE_MAX_VEC: {
#ifdef _HAVE_GMSH_
      if ( dim != 2 ) {
        MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                    "Cannot output curvatures on other entities than faces");
      }
      result = new double[9];
      pGEntity pge = EN_whatIn(pe);
      if ( GEN_type(pge) != 2 ) {
        for (int i=0; i<9; i++) result[i] = -1.;
        return result;
      }
      double u[3][2];
      F_params((pFace)pe,u);
      for (int iV=0; iV<3; iV++)
        {
          double dirMax[3], dirMin[3], curvMax, curvMin;
          GF_curvatures((pGFace)pge, u[iV],
                        dirMax, dirMin, &curvMax, &curvMin, 1.e6);
          normalizeVec(dirMax,dirMax);
          for (int iC=0; iC<3; iC++) result[3*iV+iC] = dirMax[iC] * curvMax;
        }
      return result;
#else
      MAdMsgSgl::instance().error(__LINE__,__FILE__,"Gmsh required for curvature");
#endif
    }

    case OD_CURVATURE_MIN_VEC: {
#ifdef _HAVE_GMSH_
      if ( dim != 2 ) {
        MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                    "Cannot output curvatures on other entities than faces");
      }
      result = new double[9];
      pGEntity pge = EN_whatIn(pe);
      if ( GEN_type(pge) != 2 ) {
        for (int i=0; i<9; i++) result[i] = -1.;
        return result;
      }
      double u[3][2];
      F_params((pFace)pe,u);
      for (int iV=0; iV<3; iV++)
        {
          double dirMax[3], dirMin[3], curvMax, curvMin;
          GF_curvatures((pGFace)pge, u[iV],
                        dirMax, dirMin, &curvMax, &curvMin, 1.e6);
          normalizeVec(dirMin,dirMin);
          for (int iC=0; iC<3; iC++) result[3*iV+iC] = dirMin[iC] * curvMin;
        }
      return result;
#else
      MAdMsgSgl::instance().error(__LINE__,__FILE__,"Gmsh required for curvature");
#endif
    }

    case OD_ANISO_SF_AXIS0: 
    case OD_ANISO_SF_AXIS1: 
    case OD_ANISO_SF_AXIS2:
      {
        if ( !sf ) MAdMsgSgl::instance().error(__LINE__,__FILE__,"No size field given");
      
        void * temp; 
        pPList verts;
        int k;
        pVertex pv;
        switch (dim) {
        case 3:
          result = new double[12];
          verts = R_vertices((pRegion)pe);
          temp = 0 ; k=0;
          while ( ( pv = (pVertex)PList_next(verts,&temp) ) )
            {
              pMSize size = sf->getSize(pv);
              double h, dir[3];
              if ( type == OD_ANISO_SF_AXIS0 ) h = size->direction(0,dir);
              if ( type == OD_ANISO_SF_AXIS1 ) h = size->direction(1,dir);
              if ( type == OD_ANISO_SF_AXIS2 ) h = size->direction(2,dir);
              if (size) delete size;
              for (int iC=0; iC<3; iC++) result[3*k+iC] = h * dir[iC];
              k++;
            }
          PList_delete(verts);
          return result;
        case 2:
          result = new double[9];
          verts = F_vertices( (pFace)pe, 1);
          temp = 0 ; k=0;
          while ( ( pv = (pVertex)PList_next(verts,&temp) ) )
            {
              pMSize size = sf->getSize(pv);
              double h, dir[3];
              if ( type == OD_ANISO_SF_AXIS0 ) h = size->direction(0,dir);
              if ( type == OD_ANISO_SF_AXIS1 ) h = size->direction(1,dir);
              if ( type == OD_ANISO_SF_AXIS2 ) h = size->direction(2,dir);
              if (size) delete size;
              for (int iC=0; iC<3; iC++) result[3*k+iC] = h * dir[iC];
              k++;
            }
          PList_delete(verts);
          return result;
        default:
          MAdMsgSgl::instance().error(__LINE__,__FILE__,"Dimension < 2: %d",dim);
        }
      }

    case OD_DISTANCE: {
      if ( !sf || ( sf->getType() != LOCALSFIELD ) ) {
        MAdMsgSgl::instance().error(__LINE__,__FILE__,
                                    "A local size field is required to output distance to walls");
      }
      LocalSizeField * sf_cast = (LocalSizeField *) sf;
      if ( dim == 3 ) {
        result = new double[4];
        pRegion pr = (pRegion) pe;
        for (int iV=0; iV<4; iV++) {
          result[iV] = sf_cast->getDistance(R_vertex(pr,iV));
        }
        return result;
      }
      else {
        assert ( dim == 2 );
        result = new double[3];
        pFace pf = (pFace) pe;
        for (int iV=0; iV<3; iV++) {
          result[iV] = sf_cast->getDistance(F_vertex(pf,iV));
        }
        return result;
      }
      return NULL;
    }

    default:
      MAdMsgSgl::instance().error(__LINE__,__FILE__,"Data type not handled in switch");
    }
    return 0;
  }

  // ----------------------------------------------------------------------
  void writeFaces (const pMesh m, const pSField sf, FILE *f, MAdOutputData type)
  {
    int count = 0;
    FIter fit = M_faceIter(m);
    while (pFace pf = FIter_next(fit))
      {
        double xyz[3][3];
        F_coordP1(pf, xyz);
        double* data = getData(type,(pEntity)pf,sf,count);
        if ( getOutputDataType(type) == SCALAR ) {
          fprintf (f,"ST(%g,%g,%g,%g,%g,%g,%g,%g,%g) {%g,%g,%g};\n",
                   xyz[0][0],xyz[0][1],xyz[0][2],
                   xyz[1][0],xyz[1][1],xyz[1][2],
                   xyz[2][0],xyz[2][1],xyz[2][2],
                   data[0],data[1],data[2]);
        }
        if ( getOutputDataType(type) == VECTORIAL ) {
          fprintf (f,"VT(%g,%g,%g,%g,%g,%g,%g,%g,%g) {%g,%g,%g,%g,%g,%g,%g,%g,%g};\n",
                   xyz[0][0],xyz[0][1],xyz[0][2],
                   xyz[1][0],xyz[1][1],xyz[1][2],
                   xyz[2][0],xyz[2][1],xyz[2][2],
                   data[0],data[1],data[2],
                   data[3],data[4],data[5],
                   data[6],data[7],data[8]);
        }
        delete [] data;
        count++;
      }
    FIter_delete(fit);
  }

  // ----------------------------------------------------------------------
  void writeRegions (const pMesh m, const pSField sf, FILE *f, MAdOutputData type)
  {
    int count = 0;
    RIter rit = M_regionIter(m);
    while (pRegion pr = RIter_next(rit))
      {
        double xyz[4][3];
        R_coordP1(pr, xyz);
        double* data = getData(type,(pEntity)pr,sf,count);
        if ( getOutputDataType(type) == SCALAR ) {
          fprintf (f,"SS(%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g) {%g,%g,%g,%g};\n",
                   xyz[0][0],xyz[0][1],xyz[0][2],
                   xyz[1][0],xyz[1][1],xyz[1][2],
                   xyz[2][0],xyz[2][1],xyz[2][2],
                   xyz[3][0],xyz[3][1],xyz[3][2],
                   data[0],data[1],data[2],data[3]);
        }
        if ( getOutputDataType(type) == VECTORIAL ) {
          fprintf (f,"VS(%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g) {%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g};\n",
                   xyz[0][0],xyz[0][1],xyz[0][2],
                   xyz[1][0],xyz[1][1],xyz[1][2],
                   xyz[2][0],xyz[2][1],xyz[2][2],
                   xyz[3][0],xyz[3][1],xyz[3][2],
                   data[0],data[1],data[2],data[3],
                   data[4],data[5],data[6],data[7],
                   data[8],data[9],data[10],data[11]);
        }
        delete [] data;
        count++;
      }
    RIter_delete(rit);
  }

  // ----------------------------------------------------------------------
  void MAdGmshOutput (const pMesh m, const pSField sf, const char *fn, MAdOutputData type)
  {
    MAdResourceManager& tm   = MAdResourceManagerSgl::instance();
    double t0 = tm.getTime();

    if ( sf ) MAdMsgSgl::instance().info(-1,__FILE__,
                                         "Generating output \'%s\' on file \'%s\' with size field \'%s\'",
                                         getOutputName(type).c_str(), fn, sf->getName().c_str());
    else      MAdMsgSgl::instance().info(-1,__FILE__,
                                         "Generating output \'%s\' on file \'%s\'",
                                         getOutputName(type).c_str(), fn);
  
    FILE *f = fopen (fn, "w");
    if ( !f ) {
      cerr << "Error: could not open file " << fn << endl; throw;
    }

    fprintf (f,"View\" mesh \" {\n");

    if (M_dim(m)==2) writeFaces(m,sf,f,type);
    else             writeRegions(m,sf,f,type);

    fprintf (f,"};\n");
    fclose (f);

    MAdMsgSgl::instance().info(-1,__FILE__,"Output generated in %f seconds",
                               tm.getTime()-t0);
  }

  // ----------------------------------------------------------------------
  void MAdAttachedNodalDataOutput(const pMesh m, const char *fn, pMeshDataId id)
  {
    FILE *f = fopen (fn, "w");
    if ( !f ) {
      cerr << "Error: could not open file " << fn << endl; throw;
    }

    fprintf (f,"View\" mesh \" {\n");

    if (M_dim(m)==2) {
    
      FIter fit = M_faceIter(m);
      while (pFace pf = FIter_next(fit))
        {
          // get the coordinates
          double xyz[3][3];
          F_coordP1(pf, xyz);
        
          // get the data at nodes
          double data[3];
          pPList verts = F_vertices(pf,1);
          void* tmp = 0;
          int i = 0;
          while ( pEntity pv = PList_next(verts,&tmp) ) {
            if ( EN_getDataDbl((pEntity)pv,id,&(data[i])) ) {}
            else data[i] = -10.;
            i++;
          }
          PList_delete(verts);

          // write an element
          fprintf (f,"ST(%g,%g,%g,%g,%g,%g,%g,%g,%g) {%g,%g,%g};\n",
                   xyz[0][0],xyz[0][1],xyz[0][2],
                   xyz[1][0],xyz[1][1],xyz[1][2],
                   xyz[2][0],xyz[2][1],xyz[2][2],
                   data[0],data[1],data[2]);
        }
      FIter_delete(fit);
    }
    else {

      RIter rit = M_regionIter(m);
      while (pRegion pr = RIter_next(rit))
        {
          // get the coordinates
          double xyz[4][3];
          R_coordP1(pr, xyz);
        
          // get the data at nodes
          double data[4];
          pPList verts = R_vertices(pr);
          void* tmp = 0;
          int i = 0;
          while ( pEntity pv = PList_next(verts,&tmp) ) {
            if ( EN_getDataDbl(pv,id,&(data[i])) ) {}
            else data[i] = -10.;
            i++;
          }
          PList_delete(verts);

          // write an element
          fprintf (f,"SS(%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g) {%g,%g,%g,%g};\n",
                   xyz[0][0],xyz[0][1],xyz[0][2],
                   xyz[1][0],xyz[1][1],xyz[1][2],
                   xyz[2][0],xyz[2][1],xyz[2][2],
                   xyz[3][0],xyz[3][1],xyz[3][2],
                   data[0],data[1],data[2],data[3]);
        }
      RIter_delete(rit);
  
    }

    fprintf (f,"};\n");
    fclose (f);
  }

  // ----------------------------------------------------------------------
  void MAdAttachedNodalDataVecOutput(const pMesh m, const char *fn, pMeshDataId id)
  {
    FILE *f = fopen (fn, "w");
    if ( !f ) {
      cerr << "Error: could not open file " << fn << endl; throw;
    }

    fprintf (f,"View\" mesh \" {\n");

    if (M_dim(m)==2) {
    
      FIter fit = M_faceIter(m);
      while (pFace pf = FIter_next(fit))
        {
          // get the coordinates
          double xyz[3][3];
          F_coordP1(pf, xyz);
        
          // get the data at nodes
          double data[3][3];
          pPList verts = F_vertices(pf,1);
          void* tmp = 0;
          int i = 0;
          while ( pEntity pv = PList_next(verts,&tmp) ) {
            void * tmpDat = NULL;
            if ( EN_getDataPtr((pEntity)pv,id,&tmpDat) )
              {
                //  vector<double> * vec = (vector<double>*) tmpDat;
                //  for (int j=0; j<3; j++)  data[i][j] = (*vec)[j];
                for (int j=0; j<3; j++)  data[i][j] = ((double*)tmpDat)[j];
              }
            else
              {
                for (int j=0; j<3; j++)  data[i][j] = -10.;
              }
            i++;
          }
          PList_delete(verts);

          // write an element
          fprintf (f,"VT(%g,%g,%g,%g,%g,%g,%g,%g,%g) {%g,%g,%g,%g,%g,%g,%g,%g,%g};\n",
                   xyz[0][0],xyz[0][1],xyz[0][2],
                   xyz[1][0],xyz[1][1],xyz[1][2],
                   xyz[2][0],xyz[2][1],xyz[2][2],
                   data[0][0],data[0][1],data[0][2],
                   data[1][0],data[1][1],data[1][2],
                   data[2][0],data[2][1],data[2][2]);
        }
      FIter_delete(fit);
    }
    else {

      RIter rit = M_regionIter(m);
      while (pRegion pr = RIter_next(rit))
        {
          // get the coordinates
          double xyz[4][3];
          R_coordP1(pr, xyz);
        
          // get the data at nodes
          double data[4][3];
          pPList verts = R_vertices(pr);
          void* tmp = 0;
          int i = 0;
          while ( pEntity pv = PList_next(verts,&tmp) ) {
            void * tmpDat = NULL;
            if ( EN_getDataPtr((pEntity)pv,id,&tmpDat) )
              {
                //  vector<double> * vec = (vector<double>*) tmpDat;
                //  for (int j=0; j<3; j++)  data[i][j] = (*vec)[j];
                for (int j=0; j<3; j++)  data[i][j] = ((double*)tmpDat)[j];
              }
            else
              {
                for (int j=0; j<3; j++)  data[i][j] = -10.;
              }
            i++;
          }
          PList_delete(verts);

          // write an element
          fprintf (f,"VS(%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g) {%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g};\n",
                   xyz[0][0],xyz[0][1],xyz[0][2],
                   xyz[1][0],xyz[1][1],xyz[1][2],
                   xyz[2][0],xyz[2][1],xyz[2][2],
                   xyz[3][0],xyz[3][1],xyz[3][2],
                   data[0][0],data[0][1],data[0][2],
                   data[1][0],data[1][1],data[1][2],
                   data[2][0],data[2][1],data[2][2],
                   data[3][0],data[3][1],data[3][2]);
        }
      RIter_delete(rit);
  
    }

    fprintf (f,"};\n");
    fclose (f);
  }

  // -------------------------------------------------------------------
  void printPosRegions(const set<pRegion> regs, string fn, MAdOutputData type, 
                       const pSField sf, int id)
  {
    FILE *f = fopen (fn.c_str(), "w");
    fprintf (f,"View\" mesh \" {\n");

    set<pRegion>::const_iterator rIter = regs.begin();
    for (; rIter != regs.end(); rIter++) {
      double* data = getData(type,(pEntity)(*rIter),sf);
      double xyz[4][3];
      R_coordP1(*rIter, xyz);
      fprintf (f,"SS(%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g) {%g,%g,%g,%g};\n",
               xyz[0][0],xyz[0][1],xyz[0][2],
               xyz[1][0],xyz[1][1],xyz[1][2],
               xyz[2][0],xyz[2][1],xyz[2][2],
               xyz[3][0],xyz[3][1],xyz[3][2],
               data[0],data[1],data[2],data[3]);
      delete [] data;
    }

    fprintf (f,"};\n");
    fclose (f);
  }

  // -------------------------------------------------------------------
  void printPosEntities(const pPList ents, string fn, MAdOutputData type, 
                        const pSField sf, int id)
  {
    FILE *f = fopen (fn.c_str(), "w");
    fprintf (f,"View\" mesh \" {\n");

    void *iter=0;
    while( pEntity ent = PList_next(ents,&iter) ) {
      double* data = getData(type,ent,sf);
      switch ( EN_type(ent) ) {
      case 0: 
        {
          double xyz[3];
          V_coord((pVertex)ent, xyz);
          fprintf (f,"SP(%g,%g,%g) {%g};\n",
                   xyz[0],xyz[1],xyz[2],
                   data[0]);
          break;
        }
      case 1: 
        {
          double xyz[2][3];
          E_coordP1((pEdge)ent, xyz);
          fprintf (f,"SL(%g,%g,%g,%g,%g,%g) {%g,%g};\n",
                   xyz[0][0],xyz[0][1],xyz[0][2],
                   xyz[1][0],xyz[1][1],xyz[1][2],
                   data[0],data[1]);
          break;
        }
      case 2: 
        {
          double xyz[3][3];
          F_coordP1((pFace)ent, xyz);
          fprintf (f,"ST(%g,%g,%g,%g,%g,%g,%g,%g,%g) {%g,%g,%g};\n",
                   xyz[0][0],xyz[0][1],xyz[0][2],
                   xyz[1][0],xyz[1][1],xyz[1][2],
                   xyz[2][0],xyz[2][1],xyz[2][2],
                   data[0],data[1],data[2]);
          break;
        }
      case 3: 
        {
          double xyz[4][3];
          R_coordP1((pRegion)ent, xyz);
          fprintf (f,"SS(%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g,%g) {%g,%g,%g,%g};\n",
                   xyz[0][0],xyz[0][1],xyz[0][2],
                   xyz[1][0],xyz[1][1],xyz[1][2],
                   xyz[2][0],xyz[2][1],xyz[2][2],
                   xyz[3][0],xyz[3][1],xyz[3][2],
                   data[0],data[1],data[2],data[3]);
          break;
        }
      }
      delete [] data;
    }

    fprintf (f,"};\n");
    fclose (f);
  }

  // -------------------------------------------------------------------

}