Salome HOME
changes for mpi compilation
[modules/smesh.git] / src / MEDWrapper / V2_2 / MED_V2_2_Wrapper.cxx
1 // Copyright (C) 2007-2014  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
10 //
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 // Lesser General Public License for more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22
23 #include "MED_V2_2_Wrapper.hxx"
24 #include "MED_Algorithm.hxx"
25 #include "MED_Utilities.hxx"
26
27 #include <med.h>
28 #include <med_err.h>
29
30 #ifdef _DEBUG_
31 static int MYDEBUG = 0;
32 // #else
33 // static int MYDEBUG = 0;
34 #endif
35
36
37
38 namespace MED
39 {
40   template<>
41   TInt
42   GetDESCLength<eV2_2>()
43   {
44     return 200;
45   }
46
47   template<>
48   TInt
49   GetIDENTLength<eV2_2>()
50   {
51     return 8;
52   }
53
54   template<>
55   TInt
56   GetNOMLength<eV2_2>()
57   {
58     return 64;
59   }
60
61   template<>
62   TInt
63   GetLNOMLength<eV2_2>()
64   {
65     return 80;
66   }
67
68   template<>
69   TInt
70   GetPNOMLength<eV2_2>()
71   {
72     return 16;
73   }
74
75   template<>
76   void
77   GetVersionRelease<eV2_2>(TInt& majeur, TInt& mineur, TInt& release)
78   {
79     majeur=MED_MAJOR_NUM;
80     mineur=MED_MINOR_NUM;
81     release=MED_RELEASE_NUM;
82   }
83
84   template<>
85   TInt
86   GetNbConn<eV2_2>(EGeometrieElement typmai,
87                    EEntiteMaillage typent,
88                    TInt mdim)
89   {
90     return typmai%100;
91   }
92
93   namespace V2_2
94   {
95
96     //---------------------------------------------------------------
97     class TFile{
98       TFile();
99       TFile(const TFile&);
100       
101     public:
102       TFile(const std::string& theFileName): 
103         myCount(0),
104         myFid(0), 
105         myFileName(theFileName)
106       {}
107       
108       ~TFile()
109       { 
110         Close();
111       }
112       
113       void
114       Open(EModeAcces theMode, TErr* theErr = NULL)
115       {
116         if(myCount++ == 0){
117           const char* aFileName = myFileName.c_str();
118           myFid = MEDfileOpen(aFileName,med_access_mode(theMode));
119         }
120         if(theErr)
121           *theErr = TErr(myFid);
122         else if(myFid < 0)
123           EXCEPTION(std::runtime_error,"TFile - MEDfileOpen('"<<myFileName<<"',"<<theMode<<")");
124       }
125
126       const TIdt& Id() const 
127       { 
128         if(myFid < 0)
129           EXCEPTION(std::runtime_error,"TFile - GetFid() < 0");
130         return myFid;
131       }
132
133       void Close()
134       { 
135         if(--myCount == 0)
136           MEDfileClose(myFid);
137       }
138
139     protected:
140       TInt myCount;
141       TIdt myFid;
142       std::string myFileName;
143     };
144
145
146     //---------------------------------------------------------------
147     class TFileWrapper
148     {
149       PFile myFile;
150
151     public:
152       TFileWrapper(const PFile& theFile, EModeAcces theMode, TErr* theErr = NULL): 
153         myFile(theFile)
154       {
155         myFile->Open(theMode,theErr);
156       }
157       
158       ~TFileWrapper()
159       {
160         myFile->Close();
161       }
162     };
163
164
165     //---------------------------------------------------------------
166     TVWrapper::TVWrapper(const std::string& theFileName): 
167       myFile(new TFile(theFileName))
168     {}
169     
170     
171     //----------------------------------------------------------------------------
172     TInt
173     TVWrapper
174     ::GetNbMeshes(TErr* theErr)
175     {
176       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
177       
178       if(theErr && *theErr < 0)
179         return -1;
180       
181       return MEDnMesh(myFile->Id());
182     }
183     
184     
185     //----------------------------------------------------------------------------
186     void
187     TVWrapper
188     ::GetMeshInfo(TInt theMeshId, 
189                   MED::TMeshInfo& theInfo,
190                   TErr* theErr)
191     {
192       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
193       
194       if(theErr && *theErr < 0)
195         return;
196       
197       TValueHolder<TString, char> aMeshName(theInfo.myName);
198       TValueHolder<TInt, med_int> aDim(theInfo.myDim);
199       TValueHolder<TInt, med_int> aSpaceDim(theInfo.mySpaceDim);
200       TValueHolder<EMaillage, med_mesh_type> aType(theInfo.myType);
201       char dtunit[MED_SNAME_SIZE+1];
202       med_sorting_type sorttype;
203       med_int nstep;
204       med_axis_type at;
205       int naxis=MEDmeshnAxis(myFile->Id(),theMeshId);
206       char *axisname=new char[naxis*MED_SNAME_SIZE+1];
207       char *axisunit=new char[naxis*MED_SNAME_SIZE+1];
208       TErr aRet = MEDmeshInfo(myFile->Id(),
209                               theMeshId,
210                               &aMeshName,
211                               &aSpaceDim,
212                               &aDim,
213                               &aType,
214                               &theInfo.myDesc[0],
215                               dtunit,
216                               &sorttype,
217                               &nstep,
218                               &at,
219                               axisname,
220                               axisunit);
221       delete [] axisname;
222       delete [] axisunit;
223       if(aRet < 0)
224         EXCEPTION(std::runtime_error,"GetMeshInfo - MEDmeshInfo(...)");
225     }
226     
227     
228     //----------------------------------------------------------------------------
229     void
230     TVWrapper
231     ::SetMeshInfo(const MED::TMeshInfo& theInfo,
232                   EModeAcces theMode,
233                   TErr* theErr)
234     {
235       TFileWrapper aFileWrapper(myFile,theMode,theErr);
236       
237       if(theErr && *theErr < 0)
238         return;
239       
240       MED::TMeshInfo& anInfo = const_cast<MED::TMeshInfo&>(theInfo);
241       
242       TValueHolder<TString, char> aMeshName(anInfo.myName);
243       TValueHolder<TInt, med_int> aDim(anInfo.myDim);
244       TValueHolder<TInt, med_int> aSpaceDim(anInfo.mySpaceDim);
245       TValueHolder<EMaillage, med_mesh_type> aType(anInfo.myType);
246       TValueHolder<TString, char> aDesc(anInfo.myDesc);
247
248       char *nam=new char[aSpaceDim*MED_SNAME_SIZE+1];
249       std::fill(nam,nam+aSpaceDim*MED_SNAME_SIZE+1,'\0');
250       char *unit=new char[aSpaceDim*MED_SNAME_SIZE+1];
251       std::fill(unit,unit+aSpaceDim*MED_SNAME_SIZE+1,'\0');
252       TErr aRet = MEDmeshCr(myFile->Id(),
253                             &aMeshName,
254                             aSpaceDim,
255                             aDim,
256                             aType,
257                             &aDesc,
258                             "",
259                             MED_SORT_DTIT,
260                             MED_CARTESIAN,
261                             nam,
262                             unit);
263       delete [] nam;
264       delete [] unit;
265       
266       //if(aRet == 0)
267       //  aRet = MEDunvCr(myFile->Id(),&aMeshName);
268       
269       INITMSG(MYDEBUG,"TVWrapper::SetMeshInfo - MED_MODE_ACCES = "<<theMode<<"; aRet = "<<aRet<<std::endl);
270       
271       if(theErr) 
272         *theErr = aRet;
273       else if(aRet < 0)
274         EXCEPTION(std::runtime_error,"SetMeshInfo - MEDmeshCr(...)");
275     }
276     
277     
278     //----------------------------------------------------------------------------
279     void 
280     TVWrapper
281     ::SetMeshInfo(const MED::TMeshInfo& theInfo,
282                   TErr* theErr)
283     {
284       TErr aRet;
285       SetMeshInfo(theInfo,eLECTURE_ECRITURE,&aRet);
286       
287       if(aRet < 0)
288         SetMeshInfo(theInfo,eLECTURE_AJOUT,&aRet);
289
290       if(aRet < 0)
291         SetMeshInfo(theInfo,eCREATION,&aRet);
292
293       if(theErr)
294         *theErr = aRet;
295     }
296     
297     
298     //----------------------------------------------------------------------------
299     TInt
300     TVWrapper
301     ::GetNbFamilies(const MED::TMeshInfo& theInfo,
302                     TErr* theErr)
303     {
304       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
305       
306       if(theErr && *theErr < 0)
307         return -1;
308       
309       MED::TMeshInfo& anInfo = const_cast<MED::TMeshInfo&>(theInfo);
310       TValueHolder<TString, char> aName(anInfo.myName);
311       return MEDnFamily(myFile->Id(),&aName);
312     }
313     
314     
315     //----------------------------------------------------------------------------
316     TInt
317     TVWrapper
318     ::GetNbFamAttr(TInt theFamId, 
319                    const MED::TMeshInfo& theInfo,
320                    TErr* theErr)
321     {
322       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
323       
324       if(theErr && *theErr < 0)
325         return -1;
326       
327       MED::TMeshInfo& anInfo = const_cast<MED::TMeshInfo&>(theInfo);
328
329       TValueHolder<TString, char> aName(anInfo.myName);
330
331       return MEDnFamily23Attribute(myFile->Id(),&aName,theFamId);
332     }
333     
334     
335     //----------------------------------------------------------------------------
336     TInt
337     TVWrapper
338     ::GetNbFamGroup(TInt theFamId, 
339                     const MED::TMeshInfo& theInfo,
340                     TErr* theErr)
341     {
342       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
343       
344       if(theErr && *theErr < 0)
345         return -1;
346       
347       MED::TMeshInfo& anInfo = const_cast<MED::TMeshInfo&>(theInfo);
348
349       TValueHolder<TString, char> aName(anInfo.myName);
350
351       return MEDnFamilyGroup(myFile->Id(),&aName,theFamId);
352     }
353     
354     
355     //----------------------------------------------------------------------------
356     void
357     TVWrapper
358     ::GetFamilyInfo(TInt theFamId, 
359                     MED::TFamilyInfo& theInfo,
360                     TErr* theErr)
361     {
362       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
363       
364       if(theErr && *theErr < 0)
365         return;
366       
367       MED::TMeshInfo& aMeshInfo = *theInfo.myMeshInfo;
368       
369       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
370       TValueHolder<TString, char> aFamilyName(theInfo.myName);
371       TValueHolder<TInt, med_int> aFamilyId(theInfo.myId);
372       TValueHolder<TFamAttr, med_int> anAttrId(theInfo.myAttrId);
373       TValueHolder<TFamAttr, med_int> anAttrVal(theInfo.myAttrVal);
374       TValueHolder<TString, char> anAttrDesc(theInfo.myAttrDesc);
375       TValueHolder<TString, char> aGroupNames(theInfo.myGroupNames);
376       
377       TErr aRet = MEDfamily23Info(myFile->Id(),
378                                   &aMeshName,
379                                   theFamId,
380                                   &aFamilyName,
381                                   &anAttrId,
382                                   &anAttrVal,
383                                   &anAttrDesc,
384                                   &aFamilyId,
385                                   &aGroupNames);
386
387       if(theErr) 
388         *theErr = aRet;
389       else if(aRet < 0)
390         EXCEPTION(std::runtime_error,"GetFamilyInfo - MEDfamily23Info(...) - "<<
391                   " aMeshInfo.myName = '"<<&aMeshName<<
392                   "'; theFamId = "<<theFamId<<
393                   "; theInfo.myNbGroup = "<<theInfo.myNbGroup<<
394                   "; theInfo.myNbAttr = "<<theInfo.myNbAttr);
395     }
396     
397     
398     //----------------------------------------------------------------------------
399     void
400     TVWrapper
401     ::SetFamilyInfo(const MED::TFamilyInfo& theInfo,
402                     EModeAcces theMode,
403                     TErr* theErr)
404     {
405       TFileWrapper aFileWrapper(myFile,theMode,theErr);
406       
407       if(theErr && *theErr < 0)
408         return;
409       
410       MED::TFamilyInfo& anInfo = const_cast<MED::TFamilyInfo&>(theInfo);
411       MED::TMeshInfo& aMeshInfo = *anInfo.myMeshInfo;
412       
413       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
414       TValueHolder<TString, char> aFamilyName(anInfo.myName);
415       TValueHolder<TInt, med_int> aFamilyId(anInfo.myId);
416       TValueHolder<TFamAttr, med_int> anAttrId(anInfo.myAttrId);
417       TValueHolder<TFamAttr, med_int> anAttrVal(anInfo.myAttrVal);
418       TValueHolder<TInt, med_int> aNbAttr(anInfo.myNbAttr);
419       TValueHolder<TString, char> anAttrDesc(anInfo.myAttrDesc);
420       TValueHolder<TInt, med_int> aNbGroup(anInfo.myNbGroup);
421       TValueHolder<TString, char> aGroupNames(anInfo.myGroupNames);
422
423       TErr aRet = MEDfamilyCr(myFile->Id(),
424                               &aMeshName, 
425                               &aFamilyName,
426                               aFamilyId,
427                               aNbGroup,
428                               &aGroupNames);                 
429
430       INITMSG(MYDEBUG,"TVWrapper::SetFamilyInfo - MED_MODE_ACCES = "<<theMode<<"; aRet = "<<aRet<<std::endl);
431       
432       if(theErr) 
433         *theErr = aRet;
434       else if(aRet < 0)
435         EXCEPTION(std::runtime_error,"SetFamilyInfo - MEDfamilyCr(...)");
436     }
437     
438     
439     //----------------------------------------------------------------------------
440     void
441     TVWrapper
442     ::SetFamilyInfo(const MED::TFamilyInfo& theInfo,
443                     TErr* theErr)
444     {
445       TErr aRet;
446       SetFamilyInfo(theInfo,eLECTURE_ECRITURE,&aRet);
447       
448       if(aRet < 0)
449         SetFamilyInfo(theInfo,eLECTURE_AJOUT,&aRet);
450
451       if(theErr)
452         *theErr = aRet;
453     }
454     
455     //----------------------------------------------------------------------------
456     void
457     TVWrapper
458     ::GetNames(TElemInfo&        theInfo,
459                TInt              theNb,
460                EEntiteMaillage   theEntity, 
461                EGeometrieElement theGeom,
462                TErr*             theErr)
463     {
464       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
465       
466       if(theErr && *theErr < 0)
467         return;
468       
469       if ( theGeom == eBALL )
470         theGeom = GetBallGeom( theInfo.myMeshInfo );
471
472       MED::TMeshInfo& aMeshInfo = *theInfo.myMeshInfo;
473
474       TValueHolder<TString, char>                        aMeshName  (aMeshInfo.myName);
475       TValueHolder<TString, char>                        anElemNames(theInfo.myElemNames);
476       TValueHolder<EEntiteMaillage, med_entity_type>     anEntity   (theEntity);
477       TValueHolder<EGeometrieElement, med_geometry_type> aGeom      (theGeom);
478       
479       TErr aRet = MEDmeshEntityNameRd(myFile->Id(),
480                                       &aMeshName,
481                                       MED_NO_DT,
482                                       MED_NO_IT,
483                                       anEntity,
484                                       aGeom,
485                                       &anElemNames);
486
487       theInfo.myIsElemNames = aRet != 0? eFAUX : eVRAI ;
488
489       if(theErr)
490         *theErr = aRet;
491     }
492
493     //----------------------------------------------------------------------------
494     void
495     TVWrapper
496     ::GetNumeration(TElemInfo&        theInfo,
497                     TInt              theNb,
498                     EEntiteMaillage   theEntity, 
499                     EGeometrieElement theGeom,
500                     TErr*             theErr)
501     {
502       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
503       
504       if(theErr && *theErr < 0)
505         return;
506       
507       if ( theGeom == eBALL )
508         theGeom = GetBallGeom( theInfo.myMeshInfo );
509
510       MED::TMeshInfo& aMeshInfo = *theInfo.myMeshInfo;
511       
512       TValueHolder<TString, char>                        aMeshName(aMeshInfo.myName);
513       TValueHolder<TElemNum, med_int>                    anElemNum(theInfo.myElemNum);
514       TValueHolder<EEntiteMaillage, med_entity_type>     anEntity (theEntity);
515       TValueHolder<EGeometrieElement, med_geometry_type> aGeom    (theGeom);
516       
517       TErr aRet = MEDmeshEntityNumberRd(myFile->Id(),
518                                         &aMeshName,
519                                         MED_NO_DT,
520                                         MED_NO_IT,
521                                         anEntity,
522                                         aGeom,
523                                         &anElemNum);
524
525       theInfo.myIsElemNum = aRet != 0? eFAUX : eVRAI;
526
527       if(theErr)
528         *theErr = aRet;
529     }
530
531     //----------------------------------------------------------------------------
532     void
533     TVWrapper
534     ::GetFamilies(TElemInfo&        theInfo,
535                   TInt              theNb,
536                   EEntiteMaillage   theEntity, 
537                   EGeometrieElement theGeom,
538                   TErr*             theErr)
539     {
540       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
541       
542       if(theErr && *theErr < 0)
543         return;
544
545       if ( theGeom == eBALL )
546         theGeom = GetBallGeom( theInfo.myMeshInfo );
547       
548       MED::TMeshInfo& aMeshInfo = *theInfo.myMeshInfo;
549
550       TValueHolder<TString, char>                        aMeshName(aMeshInfo.myName);
551       TValueHolder<TElemNum, med_int>                    aFamNum  (theInfo.myFamNum);
552       TValueHolder<EEntiteMaillage, med_entity_type>     anEntity (theEntity);
553       TValueHolder<EGeometrieElement, med_geometry_type> aGeom    (theGeom);
554       
555       TErr aRet = MEDmeshEntityFamilyNumberRd(myFile->Id(),
556                                               &aMeshName,
557                                               MED_NO_DT,
558                                               MED_NO_IT,
559                                               anEntity,
560                                               aGeom,
561                                               &aFamNum);
562
563       if(aRet < 0)
564       {
565 //        if (aRet == MED_ERR_DOESNTEXIST) // --- only valid with MED3.x files
566           {
567             int aSize = (int)theInfo.myFamNum->size();
568             theInfo.myFamNum->clear();
569             theInfo.myFamNum->resize(aSize,0);
570             aRet = 0;
571           }
572 //        else
573 //          EXCEPTION(std::runtime_error,"GetGrilleInfo - MEDmeshEntityFamilyNumberRd(...) of CELLS");
574       }
575       if(theErr)
576         *theErr = aRet;
577     }
578
579
580     //----------------------------------------------------------------------------
581     void
582     TVWrapper
583     ::SetNames(const TElemInfo& theInfo,
584                EEntiteMaillage theEntity, 
585                EGeometrieElement theGeom,
586                TErr* theErr)
587     { 
588       SetNames(theInfo,eLECTURE_ECRITURE,theEntity,theGeom,theErr);
589     }
590
591
592     //----------------------------------------------------------------------------
593     void
594     TVWrapper
595     ::SetNames(const TElemInfo&  theInfo,
596                EModeAcces        theMode,
597                EEntiteMaillage   theEntity, 
598                EGeometrieElement theGeom,
599                TErr*             theErr)
600     {
601       TFileWrapper aFileWrapper(myFile,theMode,theErr);
602       
603       if(theErr && *theErr < 0)
604         return;
605
606       if ( theGeom == eBALL )
607         theGeom = GetBallGeom( theInfo.myMeshInfo );
608
609       MED::TElemInfo& anInfo = const_cast<MED::TElemInfo&>(theInfo);
610       MED::TMeshInfo& aMeshInfo = *anInfo.myMeshInfo;
611
612       TErr aRet = 0;
613       if(theInfo.myIsElemNames)
614       {
615         TValueHolder<TString, char>                        aMeshName  (aMeshInfo.myName);
616         TValueHolder<TString, char>                        anElemNames(anInfo.myElemNames);
617         TValueHolder<EEntiteMaillage, med_entity_type>     anEntity   (theEntity);
618         TValueHolder<EGeometrieElement, med_geometry_type> aGeom      (theGeom);
619       
620         aRet  = MEDmeshEntityNameWr(myFile->Id(),
621                                     &aMeshName,
622                                     MED_NO_DT,
623                                     MED_NO_IT,
624                                     anEntity,
625                                     aGeom, 
626                                     (TInt)anInfo.myElemNames->size(),
627                                     &anElemNames);
628         if(theErr) 
629           *theErr = aRet;
630         else if(aRet < 0)
631           EXCEPTION(std::runtime_error,"SetNames - MEDmeshEntityNameWr(...)");
632       }
633     }
634
635
636     //----------------------------------------------------------------------------
637     void
638     TVWrapper
639     ::SetNumeration(const TElemInfo& theInfo,
640                     EEntiteMaillage theEntity, 
641                     EGeometrieElement theGeom,
642                     TErr* theErr)
643     { 
644       SetNumeration(theInfo,eLECTURE_ECRITURE,theEntity,theGeom,theErr);
645     }
646
647
648     //----------------------------------------------------------------------------
649     void 
650     TVWrapper
651     ::SetNumeration(const TElemInfo&  theInfo,
652                     EModeAcces        theMode,
653                     EEntiteMaillage   theEntity, 
654                     EGeometrieElement theGeom,
655                     TErr*             theErr)
656     {
657       TFileWrapper aFileWrapper(myFile,theMode,theErr);
658       
659       if(theErr && *theErr < 0)
660         return;
661
662       if ( theGeom == eBALL )
663         theGeom = GetBallGeom( theInfo.myMeshInfo );
664
665       MED::TElemInfo& anInfo = const_cast<MED::TElemInfo&>(theInfo);
666       MED::TMeshInfo& aMeshInfo = *anInfo.myMeshInfo;
667
668       TErr aRet = 0;
669       if(theInfo.myIsElemNum)
670       {
671         TValueHolder<TString, char>                        aMeshName(aMeshInfo.myName);
672         TValueHolder<TElemNum, med_int>                    anElemNum(anInfo.myElemNum);
673         TValueHolder<EEntiteMaillage, med_entity_type>     anEntity (theEntity);
674         TValueHolder<EGeometrieElement, med_geometry_type> aGeom    (theGeom);
675       
676         aRet  = MEDmeshEntityNumberWr(myFile->Id(),
677                                       &aMeshName,
678                                       MED_NO_DT,
679                                       MED_NO_IT,
680                                       anEntity,
681                                       aGeom,
682                                       (TInt)anInfo.myElemNum->size(),
683                                       &anElemNum);
684         if(theErr) 
685           *theErr = aRet;
686         else if(aRet < 0)
687           EXCEPTION(std::runtime_error,"SetNumeration - MEDmeshEntityNumberWr(...)");
688       }
689     }
690
691     //----------------------------------------------------------------------------
692     void
693     TVWrapper
694     ::SetFamilies(const TElemInfo&  theInfo,
695                   EEntiteMaillage   theEntity, 
696                   EGeometrieElement theGeom,
697                   TErr*             theErr)
698     { 
699       SetFamilies(theInfo,eLECTURE_ECRITURE,theEntity,theGeom,theErr);
700     }
701
702     //----------------------------------------------------------------------------
703     void 
704     TVWrapper
705     ::SetFamilies(const TElemInfo&  theInfo,
706                   EModeAcces        theMode,
707                   EEntiteMaillage   theEntity, 
708                   EGeometrieElement theGeom,
709                   TErr*             theErr)
710     {
711       TFileWrapper aFileWrapper(myFile,theMode,theErr);
712       
713       if(theErr && *theErr < 0)
714         return;
715
716       if ( theGeom == eBALL )
717         theGeom = GetBallGeom( theInfo.myMeshInfo );
718
719       MED::TElemInfo& anInfo = const_cast<MED::TElemInfo&>(theInfo);
720       MED::TMeshInfo& aMeshInfo = *anInfo.myMeshInfo;
721
722       TValueHolder<TString, char>                        aMeshName(aMeshInfo.myName);
723       TValueHolder<TElemNum, med_int>                    aFamNum  (anInfo.myFamNum);
724       TValueHolder<EEntiteMaillage, med_entity_type>     anEntity (theEntity);
725       TValueHolder<EGeometrieElement, med_geometry_type> aGeom    (theGeom);
726       
727       TErr aRet = MEDmeshEntityFamilyNumberWr(myFile->Id(),
728                                               &aMeshName,
729                                               MED_NO_DT,
730                                               MED_NO_IT,
731                                               anEntity,
732                                               aGeom,
733                                               (TInt)anInfo.myFamNum->size(),
734                                               &aFamNum);
735       
736       if(theErr) 
737         *theErr = aRet;
738       else if(aRet < 0)
739         EXCEPTION(std::runtime_error,"SetFamilies - MEDmeshEntityFamilyNumberWr(...)");
740     }
741     
742     //----------------------------------------------------------------------------
743     TInt
744     TVWrapper
745     ::GetNbNodes(const MED::TMeshInfo& theMeshInfo,
746                  ETable theTable,
747                  TErr* theErr)
748     {
749       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
750       
751       if(theErr && *theErr < 0)
752         return -1;
753       
754       MED::TMeshInfo& aMeshInfo = const_cast<MED::TMeshInfo&>(theMeshInfo);
755       
756       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
757       TValueHolder<ETable, med_data_type > aTable(theTable);
758       med_bool chgt,trsf;
759       return MEDmeshnEntity(myFile->Id(),
760                             &aMeshName,
761                             MED_NO_DT,
762                             MED_NO_IT,
763                             MED_NODE,
764                             MED_NO_GEOTYPE,
765                             aTable,
766                             MED_NO_CMODE,
767                             &chgt,
768                             &trsf);
769     }
770     
771     
772     //----------------------------------------------------------------------------
773     void
774     TVWrapper
775     ::GetNodeInfo(MED::TNodeInfo& theInfo,
776                   TErr* theErr)
777     {
778       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
779       
780       if(theErr && *theErr < 0)
781         return;
782       
783       MED::TMeshInfo& aMeshInfo = *theInfo.myMeshInfo;
784
785       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
786       TValueHolder<TInt, med_int> aDim(aMeshInfo.myDim);
787       TValueHolder<TNodeCoord, med_float> aCoord(theInfo.myCoord);
788       TValueHolder<EModeSwitch, med_switch_mode> aModeSwitch(theInfo.myModeSwitch);
789       TValueHolder<ERepere, med_axis_type> aSystem(theInfo.mySystem);
790       TValueHolder<TString, char> aCoordNames(theInfo.myCoordNames);
791       TValueHolder<TString, char> aCoordUnits(theInfo.myCoordUnits);
792       TValueHolder<TString, char> anElemNames(theInfo.myElemNames);
793       //TValueHolder<EBooleen, med_bool> anIsElemNames(theInfo.myIsElemNames);
794       TValueHolder<TElemNum, med_int> anElemNum(theInfo.myElemNum);
795       //TValueHolder<EBooleen, med_bool> anIsElemNum(theInfo.myIsElemNum);
796       TValueHolder<TElemNum, med_int> aFamNum(theInfo.myFamNum);
797       TValueHolder<TInt, med_int> aNbElem(theInfo.myNbElem);
798
799       TErr aRet = MEDmeshNodeCoordinateRd(myFile->Id(),
800                                           &aMeshName,
801                                           MED_NO_DT,
802                                           MED_NO_IT,
803                                           aModeSwitch,
804                                           &aCoord);
805
806       TErr aRet2 =MEDmeshEntityFamilyNumberRd(myFile->Id(),
807                                   &aMeshName,
808                                   MED_NO_DT,
809                                   MED_NO_IT,
810                                   MED_NODE,
811                                   MED_NO_GEOTYPE ,
812                                   &aFamNum);
813       if (aRet2  < 0)
814       {
815 //        if (aRet2 == MED_ERR_DOESNTEXIST) // --- only valid with MED3.x files
816           {
817             int mySize = (int)theInfo.myFamNum->size();
818             theInfo.myFamNum->clear();
819             theInfo.myFamNum->resize(mySize,0);
820           }
821 //        else
822 //          EXCEPTION(std::runtime_error,"GetNodeInfo - MEDmeshEntityFamilyNumberRd(...)");
823       }
824                                   
825       if ( MEDmeshEntityNameRd(myFile->Id(),
826                           &aMeshName,
827                           MED_NO_DT,
828                           MED_NO_IT,
829                           MED_NODE,
830                           MED_NO_GEOTYPE ,
831                           &anElemNames) < 0) theInfo.myIsElemNames=eFAUX;
832       
833       if ( MEDmeshEntityNumberRd(myFile->Id(),
834                             &aMeshName,
835                             MED_NO_DT,
836                             MED_NO_IT,
837                             MED_NODE,
838                             MED_NO_GEOTYPE ,
839                             &anElemNum) < 0 ) theInfo.myIsElemNum=eFAUX;
840       
841       if(theErr) 
842         *theErr = aRet;
843       else if(aRet < 0)
844         EXCEPTION(std::runtime_error,"GetNodeInfo - MEDmeshNodeCoordinateRd(...)");
845     }
846     
847     
848     //----------------------------------------------------------------------------
849     void
850     TVWrapper
851     ::SetNodeInfo(const MED::TNodeInfo& theInfo,
852                   EModeAcces theMode,
853                   TErr* theErr)
854     {
855       TFileWrapper aFileWrapper(myFile,theMode,theErr);
856       
857       if(theErr && *theErr < 0)
858         return;
859       
860       MED::TNodeInfo& anInfo = const_cast<MED::TNodeInfo&>(theInfo);
861       MED::TMeshInfo& aMeshInfo = *anInfo.myMeshInfo;
862       
863       TValueHolder<TString, char>                aMeshName    (aMeshInfo.myName);
864       TValueHolder<TNodeCoord, med_float>        aCoord       (anInfo.myCoord);
865       TValueHolder<EModeSwitch, med_switch_mode> aModeSwitch  (anInfo.myModeSwitch);
866       TValueHolder<ERepere, med_axis_type>       aSystem      (anInfo.mySystem);
867       TValueHolder<TString, char>                aCoordNames  (anInfo.myCoordNames);
868       TValueHolder<TString, char>                aCoordUnits  (anInfo.myCoordUnits);
869       TValueHolder<TString, char>                anElemNames  (anInfo.myElemNames);
870       TValueHolder<EBooleen, med_bool>           anIsElemNames(anInfo.myIsElemNames);
871       TValueHolder<TElemNum, med_int>            anElemNum    (anInfo.myElemNum);
872       TValueHolder<EBooleen, med_bool>           anIsElemNum  (anInfo.myIsElemNum);
873       TValueHolder<TElemNum, med_int>            aFamNum      (anInfo.myFamNum);
874       TValueHolder<TInt, med_int>                aNbElem      (anInfo.myNbElem);
875
876       TErr aRet = MEDmeshNodeCoordinateWr(myFile->Id(),
877                                           &aMeshName,
878                                           MED_NO_DT,
879                                           MED_NO_IT,
880                                           MED_NO_DT,
881                                           aModeSwitch,
882                                           aNbElem,
883                                           &aCoord);
884                                           
885       MEDmeshEntityFamilyNumberWr(myFile->Id(),
886                                   &aMeshName,
887                                   MED_NO_DT,
888                                   MED_NO_IT,
889                                   MED_NODE,
890                                   MED_NO_GEOTYPE,
891                                   aNbElem,
892                                   &aFamNum);
893       if(anIsElemNames)
894         MEDmeshEntityNameWr(myFile->Id(),
895                             &aMeshName,
896                             MED_NO_DT,
897                             MED_NO_IT,
898                             MED_NODE,
899                             MED_NO_GEOTYPE,
900                             aNbElem,
901                             &anElemNames);
902       if(anIsElemNum)
903         MEDmeshEntityNumberWr(myFile->Id(),
904                               &aMeshName,
905                               MED_NO_DT,
906                               MED_NO_IT,
907                               MED_NODE,
908                               MED_NO_GEOTYPE,
909                               aNbElem,
910                               &anElemNum);
911       if(theErr) 
912         *theErr = aRet;
913       else if(aRet < 0)
914         EXCEPTION(std::runtime_error,"SetNodeInfo - MEDmeshNodeCoordinateWr(...)");
915     }
916     
917     
918     //----------------------------------------------------------------------------
919     void
920     TVWrapper
921     ::SetNodeInfo(const MED::TNodeInfo& theInfo,
922                   TErr* theErr)
923     {
924       TErr aRet;
925       SetNodeInfo(theInfo,eLECTURE_ECRITURE,&aRet);
926       
927       if(aRet < 0)
928         SetNodeInfo(theInfo,eLECTURE_AJOUT,&aRet);
929
930       if(theErr) 
931         *theErr = aRet;
932     }
933     
934
935     //-----------------------------------------------------------------
936     void
937     TVWrapper
938     ::GetPolygoneInfo(MED::TPolygoneInfo& theInfo,
939                       TErr* theErr)
940     {
941       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
942
943       if(theErr && *theErr < 0)
944         return;
945
946       MED::TMeshInfo& aMeshInfo = *theInfo.myMeshInfo;
947
948       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
949       TValueHolder<TElemNum, med_int> anIndex(theInfo.myIndex);
950       TInt aNbElem = (TInt)theInfo.myElemNum->size();
951       TValueHolder<TElemNum, med_int> aConn(theInfo.myConn);
952       TValueHolder<EEntiteMaillage, med_entity_type> anEntity(theInfo.myEntity);
953       TValueHolder<EConnectivite, med_connectivity_mode> aConnMode(theInfo.myConnMode);
954
955       TErr aRet;
956       aRet = MEDmeshPolygonRd(myFile->Id(),
957                               &aMeshName,
958                               MED_NO_DT,
959                               MED_NO_IT,
960                               anEntity,
961                               aConnMode,
962                               &anIndex,
963                               &aConn);
964
965       if(theErr) 
966         *theErr = aRet;
967       else if(aRet < 0)
968         EXCEPTION(std::runtime_error,"GetPolygoneInfo - MEDmeshPolygonRd(...)");
969
970       if(theInfo.myIsElemNames){
971         GetNames(theInfo,aNbElem,theInfo.myEntity,ePOLYGONE,&aRet);
972         if(theErr) 
973           *theErr = aRet;
974       }
975
976       if(theInfo.myIsElemNum){
977         GetNumeration(theInfo,aNbElem,theInfo.myEntity,ePOLYGONE,&aRet);
978         if(theErr) 
979           *theErr = aRet;
980       }
981
982       GetFamilies(theInfo,aNbElem,theInfo.myEntity,ePOLYGONE,&aRet);
983       if(theErr) 
984         *theErr = aRet;
985     }
986     
987     //----------------------------------------------------------------------------
988     void
989     TVWrapper
990     ::SetPolygoneInfo(const MED::TPolygoneInfo& theInfo,
991                       TErr* theErr)
992     {
993       SetPolygoneInfo(theInfo,eLECTURE_ECRITURE,theErr);
994     }
995     
996     //----------------------------------------------------------------------------
997     void 
998     TVWrapper
999     ::SetPolygoneInfo(const MED::TPolygoneInfo& theInfo,
1000                       EModeAcces theMode,
1001                       TErr* theErr)
1002     {
1003       TFileWrapper aFileWrapper(myFile,theMode,theErr);
1004       
1005       if(theErr && *theErr < 0)
1006         return;
1007
1008       MED::TPolygoneInfo& anInfo = const_cast<MED::TPolygoneInfo&>(theInfo);
1009       MED::TMeshInfo& aMeshInfo = *anInfo.myMeshInfo;
1010
1011       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
1012       TValueHolder<TElemNum, med_int> anIndex(anInfo.myIndex);
1013       TValueHolder<TElemNum, med_int> aConn(anInfo.myConn);
1014       TValueHolder<EEntiteMaillage, med_entity_type> anEntity(anInfo.myEntity);
1015       TValueHolder<EConnectivite, med_connectivity_mode> aConnMode(anInfo.myConnMode);
1016
1017       TErr aRet = MEDmeshPolygonWr(myFile->Id(),
1018                                    &aMeshName,
1019                                    MED_NO_DT,
1020                                    MED_NO_IT,
1021                                    MED_UNDEF_DT,
1022                                    anEntity,
1023                                    aConnMode,
1024                                    anInfo.myNbElem + 1,
1025                                    &anIndex,
1026                                    &aConn);
1027       
1028       if(theErr) 
1029         *theErr = aRet;
1030       else if(aRet < 0)
1031         EXCEPTION(std::runtime_error,"SetPolygoneInfo - MEDmeshPolygonWr(...)");
1032       
1033       SetNames(anInfo,theInfo.myEntity,ePOLYGONE,&aRet);
1034       if(theErr) 
1035         *theErr = aRet;
1036       
1037       SetNumeration(anInfo,theInfo.myEntity,ePOLYGONE,&aRet);
1038       if(theErr) 
1039         *theErr = aRet;
1040       
1041       SetFamilies(anInfo,theInfo.myEntity,ePOLYGONE,&aRet);
1042       if(theErr) 
1043         *theErr = aRet;
1044     }
1045
1046     //----------------------------------------------------------------------------
1047     TInt 
1048     TVWrapper
1049     ::GetNbPolygones(const MED::TMeshInfo& theMeshInfo, 
1050                      EEntiteMaillage theEntity, 
1051                      EGeometrieElement theGeom, 
1052                      EConnectivite theConnMode,
1053                      TErr* theErr)
1054     {
1055       return GetNbCells(theMeshInfo,theEntity,theGeom,theConnMode,theErr);
1056     }
1057     
1058     //----------------------------------------------------------------------------
1059     TInt 
1060     TVWrapper
1061     ::GetPolygoneConnSize(const MED::TMeshInfo& theMeshInfo, 
1062                           EEntiteMaillage theEntity, 
1063                           EGeometrieElement theGeom, 
1064                           EConnectivite theConnMode,
1065                           TErr* theErr)
1066     {
1067       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1068
1069       if(theErr && *theErr < 0)
1070         return 0;
1071
1072       MED::TMeshInfo& aMeshInfo = const_cast<MED::TMeshInfo&>(theMeshInfo);
1073       
1074       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
1075       med_int aTaille = 0;
1076       med_bool chgt,trsf;
1077       aTaille=MEDmeshnEntity(myFile->Id(),
1078                              &aMeshName,
1079                              MED_NO_DT,
1080                              MED_NO_IT,
1081                              med_entity_type(theEntity),
1082                              MED_POLYGON,
1083                              MED_CONNECTIVITY,
1084                              med_connectivity_mode(theConnMode),
1085                              &chgt,
1086                              &trsf);
1087
1088       
1089       if(aTaille < 0)
1090         EXCEPTION(std::runtime_error,"GetPolygoneInfo - MEDmeshnEntity(...)");
1091
1092       return TInt(aTaille);
1093     }
1094
1095     //-----------------------------------------------------------------
1096     void 
1097     TVWrapper
1098     ::GetPolyedreInfo(TPolyedreInfo& theInfo,
1099                       TErr* theErr)
1100     {
1101       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1102
1103       if(theErr && *theErr < 0)
1104         return;
1105
1106       MED::TMeshInfo& aMeshInfo = *theInfo.myMeshInfo;
1107
1108       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
1109       TInt aNbElem = (TInt)theInfo.myElemNum->size();
1110       TValueHolder<TElemNum, med_int> anIndex(theInfo.myIndex);
1111       TValueHolder<TElemNum, med_int> aFaces(theInfo.myFaces);
1112       TValueHolder<TElemNum, med_int> aConn(theInfo.myConn);
1113       TValueHolder<EConnectivite, med_connectivity_mode> aConnMode(theInfo.myConnMode);
1114
1115       TErr aRet;
1116       aRet = MEDmeshPolyhedronRd(myFile->Id(),
1117                                  &aMeshName,
1118                                  MED_NO_DT,
1119                                  MED_NO_IT,
1120                                  MED_CELL,
1121                                  aConnMode,
1122                                  &anIndex,
1123                                  &aFaces,
1124                                  &aConn);
1125
1126       if(theErr) 
1127         *theErr = aRet;
1128       else if(aRet < 0)
1129         EXCEPTION(std::runtime_error,"GetPolygoneInfo - MEDmeshPolyhedronRd(...)");
1130
1131       if(theInfo.myIsElemNames){
1132         GetNames(theInfo,aNbElem,theInfo.myEntity,ePOLYEDRE,&aRet);
1133         if(theErr)
1134           *theErr = aRet;
1135       }
1136
1137       if(theInfo.myIsElemNum){
1138         GetNumeration(theInfo,aNbElem,theInfo.myEntity,ePOLYEDRE,&aRet);
1139         if(theErr) 
1140           *theErr = aRet;
1141       }
1142
1143       GetFamilies(theInfo,aNbElem,theInfo.myEntity,ePOLYEDRE,&aRet);
1144       if(theErr) 
1145         *theErr = aRet;
1146     }
1147
1148     //----------------------------------------------------------------------------
1149     void
1150     TVWrapper
1151     ::SetPolyedreInfo(const TPolyedreInfo& theInfo,
1152                       TErr* theErr)
1153     {
1154       SetPolyedreInfo(theInfo,eLECTURE_ECRITURE,theErr);
1155     }
1156     
1157     //----------------------------------------------------------------------------
1158     void 
1159     TVWrapper
1160     ::SetPolyedreInfo(const MED::TPolyedreInfo& theInfo,
1161                       EModeAcces theMode,
1162                       TErr* theErr)
1163     {
1164       TFileWrapper aFileWrapper(myFile,theMode,theErr);
1165       
1166       if(theErr && *theErr < 0)
1167         return;
1168
1169       MED::TPolyedreInfo& anInfo = const_cast<MED::TPolyedreInfo&>(theInfo);
1170       MED::TMeshInfo& aMeshInfo = *anInfo.myMeshInfo;
1171
1172       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
1173       TValueHolder<TElemNum, med_int> anIndex(anInfo.myIndex);
1174       TValueHolder<TElemNum, med_int> aFaces(anInfo.myFaces);
1175       TValueHolder<TElemNum, med_int> aConn(anInfo.myConn);
1176       TValueHolder<EConnectivite, med_connectivity_mode> aConnMode(anInfo.myConnMode);
1177
1178       TErr aRet;
1179       aRet = MEDmeshPolyhedronWr(myFile->Id(),
1180                                  &aMeshName,
1181                                  MED_NO_DT,
1182                                  MED_NO_IT,
1183                                  MED_UNDEF_DT,
1184                                  MED_CELL,
1185                                  aConnMode,
1186                                  anInfo.myNbElem+1,
1187                                  &anIndex,
1188                                  (TInt)anInfo.myFaces->size(),
1189                                  &aFaces,
1190                                  &aConn);
1191       
1192       if(theErr) 
1193         *theErr = aRet;
1194       else if(aRet < 0)
1195         EXCEPTION(std::runtime_error,"SetPolyedreInfo - MEDmeshPolyhedronWr(...)");
1196       
1197       TValueHolder<EEntiteMaillage, med_entity_type> anEntity(anInfo.myEntity);
1198
1199       if(theInfo.myIsElemNames){
1200         TValueHolder<TString, char> anElemNames(anInfo.myElemNames);
1201         aRet  = MEDmeshEntityNameWr(myFile->Id(),
1202                                     &aMeshName,
1203                                     MED_NO_DT,
1204                                     MED_NO_IT,
1205                                     anEntity,
1206                                     MED_POLYHEDRON,
1207                                     (TInt)anInfo.myElemNames->size(),
1208                                     &anElemNames);
1209         if(theErr) 
1210           *theErr = aRet;
1211         else if(aRet < 0)
1212           EXCEPTION(std::runtime_error,"SetPolyedreInfo - MEDmeshEntityNameWr(...)");
1213       }
1214       
1215       if(theInfo.myIsElemNum){
1216         TValueHolder<TElemNum, med_int> anElemNum(anInfo.myElemNum);
1217         aRet = MEDmeshEntityNumberWr(myFile->Id(),
1218                                      &aMeshName,
1219                                      MED_NO_DT,
1220                                      MED_NO_IT,
1221                                      anEntity,
1222                                      MED_POLYHEDRON,
1223                                      (TInt)anInfo.myElemNum->size(),
1224                                      &anElemNum);
1225         if(theErr) 
1226           *theErr = aRet;
1227         else if(aRet < 0)
1228           EXCEPTION(std::runtime_error,"SetPolyedreInfo - MEDmeshEntityNumberWr(...)");
1229       }
1230       
1231       
1232       TValueHolder<TElemNum, med_int> aFamNum(anInfo.myFamNum);
1233       aRet = MEDmeshEntityFamilyNumberWr(myFile->Id(),
1234                                          &aMeshName,
1235                                          MED_NO_DT,
1236                                          MED_NO_IT,
1237                                          anEntity,
1238                                          MED_POLYHEDRON,
1239                                          (TInt)anInfo.myFamNum->size(),
1240                                          &aFamNum);
1241       
1242       if(theErr) 
1243         *theErr = aRet;
1244       else if(aRet < 0)
1245         EXCEPTION(std::runtime_error,"SetPolyedreInfo - MEDmeshEntityFamilyNumberWr(...)");
1246     }
1247
1248     //----------------------------------------------------------------------------
1249     TInt
1250     TVWrapper
1251     ::GetNbPolyedres(const MED::TMeshInfo& theMeshInfo, 
1252                      EEntiteMaillage theEntity, 
1253                      EGeometrieElement theGeom, 
1254                      EConnectivite theConnMode,
1255                      TErr* theErr)
1256     {
1257       return GetNbCells(theMeshInfo,theEntity,theGeom,theConnMode,theErr);
1258     }
1259
1260     //----------------------------------------------------------------------------
1261     void
1262     TVWrapper    ::GetPolyedreConnSize(const TMeshInfo& theMeshInfo,
1263                           TInt& theNbFaces,
1264                           TInt& theConnSize,
1265                           EConnectivite theConnMode,
1266                           TErr* theErr)
1267     {
1268       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1269
1270       if(theErr && *theErr < 0) 
1271         EXCEPTION(std::runtime_error,"GetPolyedreConnSize - (...)");
1272
1273       MED::TMeshInfo& aMeshInfo = const_cast<MED::TMeshInfo&>(theMeshInfo);
1274       
1275       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
1276       TValueHolder<EConnectivite, med_connectivity_mode> aConnMode(theConnMode);
1277       //TValueHolder<TInt, med_int> aNbFaces(theNbFaces);
1278       //TValueHolder<TInt, med_int> aConnSize(theConnSize);
1279
1280       med_bool chgt,trsf;
1281       theNbFaces = MEDmeshnEntity(myFile->Id(),
1282                                   &aMeshName,
1283                                   MED_NO_DT,
1284                                   MED_NO_IT,
1285                                   MED_CELL,
1286                                   MED_POLYHEDRON,
1287                                   MED_INDEX_NODE,
1288                                   aConnMode,
1289                                   &chgt,
1290                                   &trsf);
1291
1292       theConnSize = MEDmeshnEntity(myFile->Id(),
1293                                   &aMeshName,
1294                                   MED_NO_DT,
1295                                   MED_NO_IT,
1296                                   MED_CELL,
1297                                   MED_POLYHEDRON,
1298                                   MED_CONNECTIVITY,
1299                                   aConnMode,
1300                                   &chgt,
1301                                   &trsf);
1302
1303       if(theNbFaces < 0 || theConnSize<0)
1304         EXCEPTION(std::runtime_error,"GetPolygoneInfo - MEDmeshnEntity(...)");
1305
1306     }
1307     
1308     //-----------------------------------------------------------------
1309     TEntityInfo
1310     TVWrapper
1311     ::GetEntityInfo(const MED::TMeshInfo& theMeshInfo,
1312                     EConnectivite theConnMode,
1313                     TErr* theErr)
1314     {
1315       TEntityInfo anInfo;
1316       
1317       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1318       
1319       if(theErr && *theErr < 0)
1320         return anInfo;
1321
1322       if(theMeshInfo.GetType() == eNON_STRUCTURE) {
1323         TInt aNbElem = GetNbNodes(theMeshInfo);
1324         if(aNbElem > 0){
1325           anInfo[eNOEUD][ePOINT1] = aNbElem;
1326           const TEntity2GeomSet& anEntity2GeomSet = GetEntity2GeomSet();
1327           TEntity2GeomSet::const_iterator anIter = anEntity2GeomSet.begin();
1328           TEntity2GeomSet::const_iterator anIterEnd = anEntity2GeomSet.end();
1329           for(; anIter != anIterEnd; anIter++){
1330             const EEntiteMaillage& anEntity = anIter->first;
1331             const TGeomSet& aGeomSet = anIter->second;
1332             TGeomSet::const_iterator anIter2 = aGeomSet.begin();
1333             TGeomSet::const_iterator anIterEnd2 = aGeomSet.end();
1334             for(; anIter2 != anIterEnd2; anIter2++){
1335               const EGeometrieElement& aGeom = *anIter2;
1336               aNbElem = GetNbCells(theMeshInfo,anEntity,aGeom,theConnMode,theErr);
1337               if(aNbElem > 0) {
1338                 if ( anEntity == eSTRUCT_ELEMENT ) {
1339                   const TInt nbStructTypes = aNbElem;
1340                   for ( TInt structType = 0; structType < nbStructTypes; ++structType ) {
1341                     // check type name to keep only "MED_BALL" structured element
1342                     TValueHolder<TString, char> aMeshName((TString&) theMeshInfo.myName );
1343                     char                        geotypename[ MED_NAME_SIZE + 1] = "";
1344                     med_geometry_type           geotype;
1345                     MEDmeshEntityInfo( myFile->Id(), &aMeshName, MED_NO_DT, MED_NO_IT,
1346                                        med_entity_type(anEntity), structType+1,
1347                                        geotypename, &geotype);
1348                     if ( strcmp( geotypename, MED_BALL_NAME ) == 0 ) {
1349                       aNbElem = GetNbCells( theMeshInfo, anEntity, EGeometrieElement(geotype),
1350                                             theConnMode, theErr);
1351                       if ( aNbElem > 0 )
1352                         anInfo[anEntity][EGeometrieElement(geotype)] = aNbElem;
1353                     }
1354                   }
1355                 }
1356                 else {
1357                   anInfo[anEntity][aGeom] = aNbElem;
1358                 }
1359               }
1360             }
1361           }
1362         }
1363       } else { // eSTRUCTURE
1364         EGrilleType aGrilleType;
1365         TInt aNbNodes = 1;
1366         TInt aNbElem  = 1;
1367         TInt aNbSub   = 0;
1368         TInt aDim = theMeshInfo.GetDim();
1369         EGeometrieElement aGeom, aSubGeom;
1370         EEntiteMaillage aSubEntity = eMAILLE;
1371
1372         GetGrilleType(theMeshInfo, aGrilleType);
1373
1374         TIntVector aStruct(aDim);
1375         if(aGrilleType == eGRILLE_STANDARD)
1376         {
1377           GetGrilleStruct(theMeshInfo, aStruct, theErr);
1378         }
1379         else
1380         { // eGRILLE_CARTESIENNE and eGRILLE_POLAIRE
1381           ETable aTable[3] = { eCOOR_IND1, eCOOR_IND2, eCOOR_IND3 };
1382           for(med_int anAxis = 0; anAxis < aDim; anAxis++)
1383             aStruct[ anAxis ] = GetNbNodes(theMeshInfo, aTable[anAxis]);
1384         }
1385         for(med_int i = 0; i < aDim; i++){
1386           aNbNodes = aNbNodes * aStruct[i];
1387           aNbElem = aNbElem * (aStruct[i] - 1);
1388         }
1389         switch(aDim){
1390         case 1:
1391           aGeom = eSEG2;
1392           break;
1393         case 2:
1394           aGeom = eQUAD4;
1395           aSubGeom = eSEG2;
1396           aSubEntity = eARETE;
1397           aNbSub =
1398             (aStruct[0]  ) * (aStruct[1]-1) +
1399             (aStruct[0]-1) * (aStruct[1]  );
1400           break;
1401         case 3:
1402           aGeom = eHEXA8;
1403           aSubGeom = eQUAD4;
1404           aSubEntity = eFACE;
1405           aNbSub =
1406             (aStruct[0]  ) * (aStruct[1]-1) * (aStruct[2]-1) +
1407             (aStruct[0]-1) * (aStruct[1]  ) * (aStruct[2]-1) +
1408             (aStruct[0]-1) * (aStruct[1]-1) * (aStruct[2]  );
1409           break;
1410         }
1411         anInfo[eNOEUD][ePOINT1] = aNbNodes;
1412         anInfo[eMAILLE][aGeom] = aNbElem;
1413         if ( aDim > 1 )
1414           anInfo[aSubEntity][aSubGeom] = aNbSub;
1415       }
1416       return anInfo;
1417     }
1418     
1419     
1420     //-----------------------------------------------------------------
1421     TInt
1422     TVWrapper
1423     ::GetNbCells(const MED::TMeshInfo& theMeshInfo, 
1424                  EEntiteMaillage theEntity, 
1425                  EGeometrieElement theGeom, 
1426                  EConnectivite theConnMode,
1427                  TErr* theErr)
1428     {
1429       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1430
1431       if(theErr && *theErr < 0)
1432         return -1;
1433
1434       MED::TMeshInfo& aMeshInfo = const_cast<MED::TMeshInfo&>(theMeshInfo);
1435       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
1436       med_bool chgt,trsf;
1437       if(theGeom!=MED::ePOLYGONE && theGeom!=MED::ePOLYEDRE && theGeom != MED::eBALL)
1438       {
1439         return MEDmeshnEntity(myFile->Id(),
1440                               &aMeshName,
1441                               MED_NO_DT,
1442                               MED_NO_IT,
1443                               med_entity_type(theEntity),
1444                               med_geometry_type(theGeom),
1445                               MED_CONNECTIVITY,
1446                               med_connectivity_mode(theConnMode),
1447                               &chgt,
1448                               &trsf);
1449       }
1450       else if(theGeom==MED::ePOLYGONE)
1451       {
1452         return MEDmeshnEntity(myFile->Id(),&aMeshName,MED_NO_DT,MED_NO_IT,med_entity_type(theEntity),
1453                               MED_POLYGON,MED_INDEX_NODE,med_connectivity_mode(theConnMode),&chgt,&trsf)-1;
1454       }
1455       else if ( theGeom==MED::ePOLYEDRE )
1456       {
1457         return MEDmeshnEntity(myFile->Id(),&aMeshName,MED_NO_DT,MED_NO_IT,med_entity_type(theEntity),
1458                               MED_POLYHEDRON,MED_INDEX_FACE,med_connectivity_mode(theConnMode),&chgt,&trsf)-1;
1459       }
1460       else if ( theGeom==MED::eBALL )
1461       {
1462         return GetNbBalls( theMeshInfo );
1463       }
1464       return 0;
1465     }
1466
1467
1468     //----------------------------------------------------------------------------
1469     void
1470     TVWrapper
1471     ::GetCellInfo(MED::TCellInfo& theInfo,
1472                   TErr* theErr)
1473     {
1474       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1475
1476       if(theErr && *theErr < 0)
1477         return;
1478       
1479       MED::TMeshInfo& aMeshInfo = *theInfo.myMeshInfo;
1480
1481       TValueHolder<TString, char>                        aMeshName    (aMeshInfo.myName);
1482       TValueHolder<TElemNum, med_int>                    aConn        (theInfo.myConn);
1483       TValueHolder<EModeSwitch, med_switch_mode>         aModeSwitch  (theInfo.myModeSwitch);
1484       TValueHolder<TString, char>                        anElemNames  (theInfo.myElemNames);
1485       TValueHolder<EBooleen, med_bool>                   anIsElemNames(theInfo.myIsElemNames);
1486       TValueHolder<TElemNum, med_int>                    anElemNum    (theInfo.myElemNum);
1487       TValueHolder<EBooleen, med_bool>                   anIsElemNum  (theInfo.myIsElemNum);
1488       TValueHolder<TElemNum, med_int>                    aFamNum      (theInfo.myFamNum);
1489       TValueHolder<EBooleen, med_bool>                   anIsFamNum   (theInfo.myIsFamNum);
1490       TValueHolder<EEntiteMaillage, med_entity_type>     anEntity     (theInfo.myEntity);
1491       TValueHolder<EGeometrieElement, med_geometry_type> aGeom        (theInfo.myGeom);
1492       TValueHolder<EConnectivite, med_connectivity_mode> aConnMode    (theInfo.myConnMode);
1493
1494       TErr aRet;
1495       aRet = MEDmeshElementRd(myFile->Id(),
1496                               &aMeshName,
1497                               MED_NO_DT,
1498                               MED_NO_IT,
1499                               anEntity,
1500                               aGeom,
1501                               aConnMode,
1502                               aModeSwitch,
1503                               &aConn,
1504                               &anIsElemNames,
1505                               &anElemNames,
1506                               &anIsElemNum,
1507                               &anElemNum,
1508                               &anIsFamNum,
1509                               &aFamNum);
1510
1511       if(theErr) 
1512         *theErr = aRet;
1513       else if(aRet < 0)
1514         EXCEPTION(std::runtime_error,"GetCellInfo - MEDmeshElementRd(...)");
1515       
1516       if (anIsFamNum == MED_FALSE)
1517         {
1518           int mySize = (int) theInfo.myFamNum->size();
1519           theInfo.myFamNum->clear();
1520           theInfo.myFamNum->resize(mySize, 0);
1521         }
1522       
1523     }
1524     
1525     
1526     //----------------------------------------------------------------------------
1527     void
1528     TVWrapper
1529     ::SetCellInfo(const MED::TCellInfo& theInfo,
1530                   EModeAcces theMode,
1531                   TErr* theErr)
1532     {
1533       TFileWrapper aFileWrapper(myFile,theMode,theErr);
1534       
1535       if(theErr && *theErr < 0)
1536         return;
1537
1538       MED::TCellInfo& anInfo    = const_cast<MED::TCellInfo&>(theInfo);
1539       MED::TMeshInfo& aMeshInfo = *anInfo.myMeshInfo;
1540
1541       TValueHolder<TString, char>                        aMeshName    (aMeshInfo.myName);
1542       TValueHolder<TElemNum, med_int>                    aConn        (anInfo.myConn);
1543       TValueHolder<EModeSwitch, med_switch_mode>         aModeSwitch  (anInfo.myModeSwitch);
1544       TValueHolder<TString, char>                        anElemNames  (anInfo.myElemNames);
1545       TValueHolder<EBooleen, med_bool>                   anIsElemNames(anInfo.myIsElemNames);
1546       TValueHolder<TElemNum, med_int>                    anElemNum    (anInfo.myElemNum);
1547       TValueHolder<EBooleen, med_bool>                   anIsElemNum  (anInfo.myIsElemNum);
1548       TValueHolder<TElemNum, med_int>                    aFamNum      (anInfo.myFamNum);
1549       TValueHolder<EBooleen, med_bool>                   anIsFamNum   (anInfo.myIsFamNum);
1550       TValueHolder<EEntiteMaillage, med_entity_type>     anEntity     (anInfo.myEntity);
1551       TValueHolder<EGeometrieElement, med_geometry_type> aGeom        (anInfo.myGeom);
1552       TValueHolder<EConnectivite, med_connectivity_mode> aConnMode    (anInfo.myConnMode);
1553       TValueHolder<TInt, med_int>                        aNbElem      (anInfo.myNbElem);
1554
1555       TErr aRet;
1556       aRet = MEDmeshElementConnectivityWr(myFile->Id(),
1557                                           &aMeshName,
1558                                           MED_NO_DT,
1559                                           MED_NO_IT,
1560                                           MED_UNDEF_DT,
1561                                           anEntity,
1562                                           aGeom,
1563                                           aConnMode,
1564                                           aModeSwitch,
1565                                           aNbElem,
1566                                           &aConn);
1567
1568       MEDmeshEntityFamilyNumberWr(myFile->Id(),
1569                                   &aMeshName,
1570                                   MED_NO_DT,
1571                                   MED_NO_IT,
1572                                   anEntity,
1573                                   aGeom,
1574                                   aNbElem,
1575                                   &aFamNum);
1576       if(anIsElemNames)
1577         MEDmeshEntityNameWr(myFile->Id(),
1578                             &aMeshName,
1579                             MED_NO_DT,
1580                             MED_NO_IT,
1581                             anEntity,
1582                             aGeom,
1583                             aNbElem,
1584                             &anElemNames);
1585       if(anIsElemNum)
1586         MEDmeshEntityNumberWr(myFile->Id(),
1587                               &aMeshName,
1588                               MED_NO_DT,
1589                               MED_NO_IT,
1590                               anEntity,
1591                               aGeom,
1592                               aNbElem,
1593                               &anElemNum);
1594       if(theErr) 
1595         *theErr = aRet;
1596       else if(aRet < 0)
1597         EXCEPTION(std::runtime_error,"SetCellInfo - MEDmeshElementWr(...)");
1598     }
1599     
1600
1601     //----------------------------------------------------------------------------
1602     void
1603     TVWrapper
1604     ::SetCellInfo(const MED::TCellInfo& theInfo,
1605                   TErr* theErr)
1606     {
1607       SetCellInfo(theInfo,eLECTURE_ECRITURE,theErr);
1608     }
1609
1610     //----------------------------------------------------------------------------
1611     //! Read geom type of MED_BALL structural element
1612     EGeometrieElement TVWrapper::GetBallGeom(const TMeshInfo& theMeshInfo)
1613     {
1614       TFileWrapper aFileWrapper(myFile,eLECTURE);
1615
1616       // read med_geometry_type of "MED_BALL" element
1617       char geotypename[ MED_NAME_SIZE + 1] = MED_BALL_NAME;
1618       return EGeometrieElement( MEDstructElementGeotype( myFile->Id(), geotypename ) );
1619     }
1620
1621     //----------------------------------------------------------------------------
1622     //! Read number of balls in the Mesh
1623     TInt TVWrapper::GetNbBalls(const TMeshInfo& theMeshInfo)
1624     {
1625       TFileWrapper aFileWrapper(myFile,eLECTURE);
1626
1627       EGeometrieElement ballType = GetBallGeom( theMeshInfo );
1628       if ( ballType < 0 )
1629         return 0;
1630
1631       return GetNbCells( theMeshInfo, eSTRUCT_ELEMENT, ballType, eNOD );
1632     }
1633
1634     //----------------------------------------------------------------------------
1635     //! Read a MEDWrapped representation of MED_BALL from the MED file
1636     void TVWrapper::GetBallInfo(TBallInfo& theInfo, TErr* theErr)
1637     {
1638       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1639
1640       // check geometry of MED_BALL
1641       if ( theInfo.myGeom == eBALL )
1642       {
1643         theInfo.myGeom = GetBallGeom( *theInfo.myMeshInfo );
1644         if ( theInfo.myGeom < 0 ) {
1645           if ( !theErr )
1646             EXCEPTION(std::runtime_error,"GetBallInfo - no balls in the mesh");
1647           *theErr = theInfo.myGeom;
1648           return;
1649         }
1650       }
1651
1652       // read nodes ids
1653       GetCellInfo( theInfo );
1654
1655       // read diameters
1656       TValueHolder<TString, char>                        aMeshName (theInfo.myMeshInfo->myName);
1657       TValueHolder<EGeometrieElement, med_geometry_type> aGeom     (theInfo.myGeom);
1658       TValueHolder<TFloatVector, void>                   aDiam     (theInfo.myDiameters);
1659       char varattname[ MED_NAME_SIZE + 1] = MED_BALL_DIAMETER;
1660
1661       TErr aRet = MEDmeshStructElementVarAttRd( myFile->Id(), &aMeshName,
1662                                                 MED_NO_DT, MED_NO_IT,
1663                                                 aGeom,
1664                                                 varattname,
1665                                                 &aDiam);
1666       if ( theErr )
1667         *theErr = aRet;
1668       else if ( aRet < 0 )
1669         EXCEPTION(std::runtime_error,"GetBallInfo - pb at reading diameters");
1670     }
1671
1672
1673     //----------------------------------------------------------------------------
1674     //! Write a MEDWrapped representation of MED_BALL to the MED file
1675     void  TVWrapper::SetBallInfo(const TBallInfo& theInfo, EModeAcces theMode, TErr* theErr)
1676     {
1677       TFileWrapper aFileWrapper(myFile,theMode,theErr);
1678
1679       TErr ret;
1680       char ballsupportname[MED_NAME_SIZE+1]="BALL_SUPPORT_MESH";
1681       EGeometrieElement ballGeom = GetBallGeom( *theInfo.myMeshInfo );
1682       if ( ballGeom < 0 )
1683       {
1684         // no ball model in the file, create support mesh for it
1685         char dummyname [MED_NAME_SIZE*3+1]="";
1686         if (( ret = MEDsupportMeshCr( myFile->Id(),
1687                                       ballsupportname,
1688                                       theInfo.myMeshInfo->GetSpaceDim(),
1689                                       theInfo.myMeshInfo->GetDim(),
1690                                       "Support mesh for a ball model",
1691                                       MED_CARTESIAN,
1692                                       /*axisname=*/dummyname, /*unitname=*/dummyname)) < 0) {
1693           if ( !theErr )
1694             EXCEPTION(std::runtime_error,"SetBallInfo - MEDsupportMeshCr");
1695           *theErr = ret;
1696           return;
1697         }
1698         // write coordinates of 1 node
1699         med_float coord[3] = {0,0,0};
1700         if ((ret = MEDmeshNodeCoordinateWr(myFile->Id(),
1701                                            ballsupportname, MED_NO_DT, MED_NO_IT, 0.0,
1702                                            MED_FULL_INTERLACE, /*nnode=*/1, coord)) < 0) {
1703           if ( !theErr )
1704             EXCEPTION(std::runtime_error,"SetBallInfo - MEDmeshNodeCoordinateWr");
1705           *theErr = ret;
1706           return;
1707         }
1708         // ball model creation
1709         char geotypename[ MED_NAME_SIZE + 1] = MED_BALL_NAME;
1710         if (( ballGeom = (EGeometrieElement) MEDstructElementCr(myFile->Id(),
1711                                                                 geotypename,
1712                                                                 theInfo.myMeshInfo->GetSpaceDim(),
1713                                                                 ballsupportname,
1714                                                                 MED_NODE,MED_NONE)) < 0 ) {
1715           if ( !theErr )
1716             EXCEPTION(std::runtime_error,"SetBallInfo - MEDstructElementCr");
1717           *theErr = ret;
1718           return;
1719         }
1720         // create diameter attribute
1721         if (( ret = MEDstructElementVarAttCr(myFile->Id(),
1722                                              geotypename, MED_BALL_DIAMETER,
1723                                              MED_ATT_FLOAT64, /*ncomp=*/1)) < 0) {
1724           if ( !theErr )
1725             EXCEPTION(std::runtime_error,"SetBallInfo - MEDstructElementVarAttCr");
1726           *theErr = ret;
1727           return;
1728         }
1729       } // ballGeom < 0
1730
1731       TBallInfo& aBallInfo = ((TBallInfo&) theInfo );
1732       aBallInfo.myGeom = ballGeom;
1733
1734       // write node ids
1735       SetCellInfo(theInfo,theMode,theErr);
1736       if ( theErr && theErr < 0 )
1737         return;
1738
1739       // write diameter
1740       TValueHolder<TString, char>                        aMeshName (aBallInfo.myMeshInfo->myName);
1741       TValueHolder<EGeometrieElement, med_geometry_type> aGeom     (aBallInfo.myGeom);
1742       TValueHolder<TFloatVector, void>                   aDiam     (aBallInfo.myDiameters);
1743       ret = MEDmeshStructElementVarAttWr(myFile->Id(), &aMeshName,
1744                                          MED_NO_DT, MED_NO_IT,
1745                                          aGeom, MED_BALL_DIAMETER,
1746                                          theInfo.myNbElem, &aDiam);
1747       if ( theErr )
1748         *theErr = ret;
1749       else if ( ret < 0 )
1750         EXCEPTION(std::runtime_error,"SetBallInfo - MEDmeshStructElementVarAttWr");
1751     }
1752
1753     //----------------------------------------------------------------------------
1754     //! Write a MEDWrapped representation of MED_BALL to the MED file
1755     void  TVWrapper::SetBallInfo(const TBallInfo& theInfo, TErr* theErr)
1756     {
1757       SetBallInfo( theInfo, eLECTURE_ECRITURE, theErr );
1758     }
1759
1760     //-----------------------------------------------------------------
1761     TInt
1762     TVWrapper
1763     ::GetNbFields(TErr* theErr)
1764     {
1765       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1766       
1767       if(theErr && *theErr < 0)
1768         return -1;
1769       
1770       return MEDnField(myFile->Id());
1771     }
1772     
1773     
1774     //----------------------------------------------------------------------------
1775     TInt
1776     TVWrapper
1777     ::GetNbComp(TInt theFieldId,
1778                 TErr* theErr)
1779     {
1780       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1781       
1782       if(theErr && *theErr < 0)
1783         return -1;
1784       
1785       return MEDfieldnComponent(myFile->Id(),theFieldId);
1786     }
1787     
1788     
1789     //----------------------------------------------------------------------------
1790     void
1791     TVWrapper
1792     ::GetFieldInfo(TInt theFieldId, 
1793                    MED::TFieldInfo& theInfo,
1794                    TErr* theErr)
1795     {
1796       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1797       
1798       if(theErr && *theErr < 0)
1799         return;
1800       
1801       TString aFieldName(256); // Protect from memory problems with too long names
1802       TValueHolder<ETypeChamp, med_field_type> aType(theInfo.myType);
1803       TValueHolder<TString, char> aCompNames(theInfo.myCompNames);
1804       TValueHolder<TString, char> anUnitNames(theInfo.myUnitNames);
1805       MED::TMeshInfo& aMeshInfo = theInfo.myMeshInfo;
1806       
1807       TErr aRet;
1808       med_bool local;
1809       char dtunit[MED_SNAME_SIZE+1];
1810           char local_mesh_name[MED_NAME_SIZE+1]="";
1811       med_int nbofstp;
1812       theInfo.myNbComp = MEDfieldnComponent(myFile->Id(),theFieldId);
1813       aRet = MEDfieldInfo(myFile->Id(),
1814                           theFieldId,
1815                           &aFieldName[0],
1816                           local_mesh_name,
1817                           &local,
1818                           &aType,
1819                           &aCompNames,
1820                           &anUnitNames,
1821                           dtunit,
1822                           &nbofstp);
1823
1824           if(strcmp(&aMeshInfo.myName[0],local_mesh_name) != 0 ) {
1825                   if(theErr)
1826                         *theErr = -1;
1827                   return;
1828           }
1829
1830       theInfo.SetName(aFieldName);
1831
1832       if(theErr)
1833         *theErr = aRet;
1834       else if(aRet < 0)
1835         EXCEPTION(std::runtime_error,"GetFieldInfo - MEDfieldInfo(...)");
1836     }
1837     
1838     
1839     //----------------------------------------------------------------------------
1840     void
1841     TVWrapper
1842     ::SetFieldInfo(const MED::TFieldInfo& theInfo,
1843                    EModeAcces theMode,
1844                    TErr* theErr)
1845     {
1846       TFileWrapper aFileWrapper(myFile,theMode,theErr);
1847       
1848       if(theErr && *theErr < 0)
1849         return;
1850       
1851       MED::TFieldInfo& anInfo = const_cast<MED::TFieldInfo&>(theInfo);
1852       
1853       TValueHolder<TString, char> aFieldName(anInfo.myName);
1854       TValueHolder<ETypeChamp, med_field_type> aType(anInfo.myType);
1855       TValueHolder<TString, char> aCompNames(anInfo.myCompNames);
1856       TValueHolder<TString, char> anUnitNames(anInfo.myUnitNames);
1857       MED::TMeshInfo& aMeshInfo = anInfo.myMeshInfo;
1858       TErr aRet;
1859       char dtunit[MED_SNAME_SIZE+1];
1860       std::fill(dtunit,dtunit+MED_SNAME_SIZE+1,'\0');
1861       aRet = MEDfieldCr(myFile->Id(),
1862                         &aFieldName,
1863                         aType,
1864                         anInfo.myNbComp,
1865                         &aCompNames,
1866                         &anUnitNames,
1867                         dtunit,
1868                         &aMeshInfo.myName[0]);
1869       if(theErr) 
1870         *theErr = aRet;
1871       else if(aRet < 0)
1872         EXCEPTION(std::runtime_error,"SetFieldInfo - MEDfieldCr(...)");
1873     }
1874     
1875     
1876     //----------------------------------------------------------------------------
1877     void
1878     TVWrapper
1879     ::SetFieldInfo(const MED::TFieldInfo& theInfo,
1880                    TErr* theErr)
1881     {
1882       TErr aRet;
1883       SetFieldInfo(theInfo,eLECTURE_ECRITURE,&aRet);
1884       
1885       if(aRet < 0)
1886         SetFieldInfo(theInfo,eLECTURE_AJOUT,&aRet);
1887
1888       if(theErr) 
1889         *theErr = aRet;
1890     }
1891     
1892     
1893     //----------------------------------------------------------------------------
1894     TInt
1895     TVWrapper
1896     ::GetNbGauss(TErr* theErr)
1897     {
1898       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1899       
1900       if(theErr && *theErr < 0)
1901         return -1;
1902       
1903       return MEDnLocalization(myFile->Id());
1904     }
1905
1906
1907     //----------------------------------------------------------------------------
1908     TGaussInfo::TInfo
1909     TVWrapper
1910     ::GetGaussPreInfo(TInt theId, 
1911                       TErr* theErr)
1912     {
1913       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1914
1915       if(theErr && *theErr < 0)
1916         return TGaussInfo::TInfo( TGaussInfo::TKey(ePOINT1,""),0 );
1917       
1918       med_int aNbGaussPoints = med_int();
1919       TVector<char> aName(GetNOMLength<eV2_2>()+1);
1920       med_geometry_type aGeom = MED_NONE;
1921
1922       TErr aRet;
1923       med_int dim;
1924       char geointerpname[MED_NAME_SIZE+1]="";
1925       char ipointstructmeshname[MED_NAME_SIZE+1]="";
1926       med_int nsectionmeshcell;
1927       med_geometry_type sectiongeotype;
1928       aRet = MEDlocalizationInfo (myFile->Id(),
1929                                   theId,
1930                                   &aName[0],
1931                                   &aGeom,
1932                                   &dim,
1933                                   &aNbGaussPoints,
1934                                   geointerpname,
1935                                   ipointstructmeshname,
1936                                   &nsectionmeshcell,
1937                                   &sectiongeotype);
1938       if(theErr) 
1939         *theErr = aRet;
1940       else if(aRet < 0)
1941         EXCEPTION(std::runtime_error,"GetGaussPreInfo - MEDlocalizationInfo(...)");
1942       return TGaussInfo::TInfo(TGaussInfo::TKey(EGeometrieElement(aGeom),&aName[0]),
1943                                TInt(aNbGaussPoints));
1944     }
1945
1946
1947     //----------------------------------------------------------------------------
1948     void
1949     TVWrapper
1950     ::GetGaussInfo(TInt theId, 
1951                    TGaussInfo& theInfo,
1952                    TErr* theErr)
1953     {
1954       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1955       
1956       if(theErr && *theErr < 0)
1957         return;
1958       
1959       TValueHolder<TNodeCoord, med_float> aRefCoord(theInfo.myRefCoord);
1960       TValueHolder<TNodeCoord, med_float> aGaussCoord(theInfo.myGaussCoord);
1961       TValueHolder<TWeight, med_float> aWeight(theInfo.myWeight);
1962       TValueHolder<EModeSwitch, med_switch_mode> aModeSwitch(theInfo.myModeSwitch);
1963       TValueHolder<TString, char> aGaussName(theInfo.myName);
1964
1965       TErr aRet;
1966       aRet = MEDlocalizationRd(myFile->Id(),
1967                                &aGaussName,
1968                                aModeSwitch,
1969                                &aRefCoord,
1970                                &aGaussCoord,
1971                                &aWeight);
1972
1973       if(theErr) 
1974         *theErr = aRet;
1975       else if(aRet < 0)
1976         EXCEPTION(std::runtime_error,"GetGaussInfo - MEDlocalizationRd(...)");
1977     }
1978
1979
1980     //----------------------------------------------------------------------------
1981     TInt
1982     TVWrapper
1983     ::GetNbProfiles(TErr* theErr)
1984     {
1985       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1986       
1987       if(theErr && *theErr < 0)
1988         return -1;
1989       
1990       return MEDnProfile(myFile->Id());
1991     }
1992
1993     TProfileInfo::TInfo
1994     TVWrapper
1995     ::GetProfilePreInfo(TInt theId, 
1996                         TErr* theErr)
1997     {
1998       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
1999
2000       if(theErr && *theErr < 0)
2001         return TProfileInfo::TInfo();
2002       
2003       med_int aSize = -1;
2004       TVector<char> aName(GetNOMLength<eV2_2>()+1);
2005
2006       TErr aRet;
2007       aRet = MEDprofileInfo(myFile->Id(),
2008                             theId,
2009                             &aName[0],
2010                             &aSize);
2011       if(theErr) 
2012         *theErr = aRet;
2013       else if(aRet < 0)
2014         EXCEPTION(std::runtime_error,"GetProfilePreInfo - MEDprofileInfo(...)");
2015       
2016       return TProfileInfo::TInfo(&aName[0],aSize);
2017     }
2018
2019     void
2020     TVWrapper
2021     ::GetProfileInfo(TInt theId, 
2022                      TProfileInfo& theInfo,
2023                      TErr* theErr)
2024     {
2025       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
2026       
2027       if(theErr && *theErr < 0)
2028         return;
2029       
2030       TProfileInfo& anInfo = const_cast<TProfileInfo&>(theInfo);
2031       TValueHolder<TElemNum, med_int> anElemNum(anInfo.myElemNum);
2032       TValueHolder<TString, char>     aProfileName(anInfo.myName);
2033
2034       TErr aRet;
2035       aRet = MEDprofileRd(myFile->Id(),
2036                           &aProfileName,
2037                           &anElemNum);
2038       if(theErr) 
2039         *theErr = aRet;
2040       else if(aRet < 0)
2041         EXCEPTION(std::runtime_error,"GetProfileInfo - MEDprofileRd(...)");
2042     }
2043
2044     void
2045     TVWrapper
2046     ::SetProfileInfo(const TProfileInfo& theInfo,
2047                      EModeAcces          theMode,
2048                      TErr*               theErr)
2049     {
2050       TFileWrapper aFileWrapper(myFile,theMode,theErr);
2051       
2052       if(theErr && *theErr < 0)
2053         return;
2054       
2055       TProfileInfo& anInfo = const_cast<TProfileInfo&>(theInfo);
2056       TValueHolder<TElemNum, med_int> anElemNum(anInfo.myElemNum);
2057       TValueHolder<TString, char>     aProfileName(anInfo.myName);
2058
2059       TErr aRet;
2060       aRet = MEDprofileWr(myFile->Id(),      // descripteur du fichier.
2061                           &aProfileName,        // tableau de valeurs du profil.
2062                           theInfo.GetSize(), // taille du profil.
2063                           &anElemNum);    // nom profil.
2064       if(theErr)
2065         *theErr = aRet;
2066       else if(aRet < 0)
2067         EXCEPTION(std::runtime_error,"SetProfileInfo - MEDprofileWr(...)");
2068     }
2069
2070     void
2071     TVWrapper
2072     ::SetProfileInfo(const TProfileInfo& theInfo,
2073                      TErr*               theErr)
2074     {
2075       TErr aRet;
2076       SetProfileInfo(theInfo,eLECTURE_ECRITURE,&aRet);
2077       
2078       if(aRet < 0)
2079         SetProfileInfo(theInfo,eLECTURE_AJOUT,&aRet);
2080
2081       if(aRet < 0)
2082         SetProfileInfo(theInfo,eCREATION,&aRet);
2083
2084       if(theErr)
2085         *theErr = aRet;
2086     }
2087
2088     //-----------------------------------------------------------------
2089     TInt
2090     TVWrapper
2091     ::GetNbTimeStamps(const MED::TFieldInfo& theInfo, 
2092                       const MED::TEntityInfo& theEntityInfo,
2093                       EEntiteMaillage& theEntity,
2094                       TGeom2Size& theGeom2Size,
2095                       TErr* theErr)
2096     {
2097       theEntity = EEntiteMaillage(-1);
2098       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
2099
2100       if(theErr){
2101         if(theEntityInfo.empty())
2102           *theErr = -1;
2103         if(*theErr < 0)
2104           return -1;
2105       }else if(theEntityInfo.empty()) 
2106         EXCEPTION(std::runtime_error,"GetNbTimeStamps - There is no any Entity on the Mesh");
2107       
2108       bool anIsPerformAdditionalCheck = GetNbMeshes() > 1;
2109
2110       theGeom2Size.clear();
2111       TInt aNbTimeStamps = 0;
2112       TIdt anId = myFile->Id();
2113
2114       MED::TFieldInfo& anInfo = const_cast<MED::TFieldInfo&>(theInfo);
2115       TValueHolder<TString, char> aFieldName(anInfo.myName);
2116       MED::TMeshInfo& aMeshInfo = anInfo.myMeshInfo;
2117
2118       // workaround for IPAL13676
2119       MED::TEntityInfo localEntityInfo = theEntityInfo;
2120       TEntityInfo::iterator anLocalIter = localEntityInfo.find(eMAILLE);
2121       if(anLocalIter != localEntityInfo.end()){
2122         localEntityInfo[eNOEUD_ELEMENT] = anLocalIter->second;
2123       }
2124         
2125       TEntityInfo::const_iterator anIter = localEntityInfo.begin();
2126       for(; anIter != localEntityInfo.end(); anIter++){
2127         med_entity_type anEntity = med_entity_type(anIter->first);
2128         const TGeom2Size& aGeom2Size = anIter->second;
2129         TGeom2Size::const_iterator anGeomIter = aGeom2Size.begin();
2130         for(; anGeomIter != aGeom2Size.end(); anGeomIter++){
2131           med_geometry_type aGeom = med_geometry_type(anGeomIter->first);
2132           char aMeshName[MED_NAME_SIZE+1];
2133           med_bool islocal;
2134           med_field_type ft;
2135           char dtunit[MED_SNAME_SIZE+1];
2136           med_int myNbComp = MEDfieldnComponentByName(anId,&aFieldName);
2137           char *cname=new char[myNbComp*MED_SNAME_SIZE+1];
2138           char *unitname=new char[myNbComp*MED_SNAME_SIZE+1];
2139           TInt aNbStamps;
2140           MEDfieldInfoByName(anId,
2141                              &aFieldName,
2142                              aMeshName,
2143                              &islocal,
2144                              &ft,
2145                              cname,
2146                              unitname,
2147                              dtunit,
2148                              &aNbStamps);
2149           delete [] cname;
2150           delete [] unitname;
2151           med_int nval = 0;
2152           med_int aNumDt;
2153           med_int aNumOrd;
2154           med_float aDt;
2155           if (aNbStamps > 0)
2156             {
2157               MEDfieldComputingStepInfo(anId,
2158                                         &aFieldName,
2159                                         1,
2160                                         &aNumDt,
2161                                         &aNumOrd,
2162                                         &aDt);
2163               char profilename[MED_NAME_SIZE+1];
2164               char locname[MED_NAME_SIZE+1];
2165               med_int profilsize;
2166               med_int aNbGauss;
2167
2168               // protection from crash (division by zero)
2169               // inside MEDfieldnValueWithProfile function
2170               // caused by the workaround for IPAL13676 (see above)
2171               if( anEntity == MED_NODE_ELEMENT && aGeom % 100 == 0 )
2172                 continue;
2173
2174               nval = MEDfieldnValueWithProfile(anId,
2175                                                &aFieldName,
2176                                                aNumDt,
2177                                                aNumOrd,
2178                                                anEntity,
2179                                                med_geometry_type(aGeom),
2180                                                1,
2181                                                MED_COMPACT_STMODE,
2182                                                profilename,
2183                                                &profilsize,
2184                                                locname,
2185                                                &aNbGauss);
2186             }
2187           bool anIsSatisfied =(nval > 0);
2188           if(anIsSatisfied){
2189             INITMSG(MYDEBUG,
2190                     "GetNbTimeStamps aNbTimeStamps = "<<aNbStamps<<
2191                     "; aGeom = "<<aGeom<<"; anEntity = "<<anEntity<<"\n");
2192             if(anIsPerformAdditionalCheck){
2193               anIsSatisfied = !strcmp(&aMeshName[0],&aMeshInfo.myName[0]);
2194               if(!anIsSatisfied){
2195                 INITMSG(MYDEBUG,
2196                         "GetNbTimeStamps aMeshName = '"<<&aMeshName[0]<<"' != "<<
2197                         "; aMeshInfo.myName = '"<<&aMeshInfo.myName[0]<<"'\n");
2198               }
2199             }
2200           }
2201           if(anIsSatisfied){
2202             theGeom2Size[EGeometrieElement(aGeom)] = anGeomIter->second;
2203             theEntity = EEntiteMaillage(anEntity);
2204             aNbTimeStamps = aNbStamps;
2205           }
2206         }
2207         if(!theGeom2Size.empty()) 
2208           break;
2209       }
2210       return aNbTimeStamps;
2211     }
2212     
2213     
2214     //----------------------------------------------------------------------------
2215     void
2216     TVWrapper
2217     ::GetTimeStampInfo(TInt theTimeStampId, 
2218                        MED::TTimeStampInfo& theInfo,
2219                        TErr* theErr)
2220     {
2221       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
2222       
2223       const TGeom2Size& aGeom2Size = theInfo.myGeom2Size;
2224       
2225       if(theErr){
2226         if(aGeom2Size.empty())
2227           *theErr = -1;
2228         if(*theErr < 0)
2229           return;
2230       }else if(aGeom2Size.empty())
2231         EXCEPTION(std::runtime_error,"GetTimeStampInfo - There is no any cell");
2232       
2233       MED::TFieldInfo& aFieldInfo = *theInfo.myFieldInfo;
2234       MED::TMeshInfo& aMeshInfo = *aFieldInfo.myMeshInfo;
2235       
2236       TValueHolder<TString, char> aFieldName(aFieldInfo.myName);
2237       TValueHolder<EEntiteMaillage, med_entity_type> anEntity(theInfo.myEntity);
2238       TValueHolder<TInt, med_int> aNumDt(theInfo.myNumDt);
2239       TValueHolder<TInt, med_int> aNumOrd(theInfo.myNumOrd);
2240       TValueHolder<TString, char> anUnitDt(theInfo.myUnitDt);
2241       TValueHolder<TFloat, med_float> aDt(theInfo.myDt);
2242       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
2243       TValueHolder<EBooleen, med_bool> anIsLocal(aFieldInfo.myIsLocal);
2244       TValueHolder<TInt, med_int> aNbRef(aFieldInfo.myNbRef);
2245
2246       TGeom2NbGauss& aGeom2NbGauss = theInfo.myGeom2NbGauss;
2247
2248       // just to get a time stamp unit (anUnitDt)
2249       med_field_type aFieldType;
2250       med_int aNbComp = MEDfieldnComponentByName(myFile->Id(), &aFieldName);
2251       char *aCompName = new char[aNbComp*MED_SNAME_SIZE+1];
2252       char *aCompUnit = new char[aNbComp*MED_SNAME_SIZE+1];
2253       TInt aNbStamps;
2254       MEDfieldInfoByName(myFile->Id(),
2255                          &aFieldName,
2256                          &aMeshName,
2257                          &anIsLocal,
2258                          &aFieldType,
2259                          aCompName,
2260                          aCompUnit,
2261                          &anUnitDt,
2262                          &aNbStamps);
2263       delete [] aCompName;
2264       delete [] aCompUnit;
2265
2266       TGeom2Size::const_iterator anIter = aGeom2Size.begin();
2267       for(; anIter != aGeom2Size.end(); anIter++){
2268         const EGeometrieElement& aGeom = anIter->first;
2269         med_int aNbGauss = -1;
2270
2271         TErr aRet;
2272         aRet = MEDfieldComputingStepInfo(myFile->Id(),
2273                                          &aFieldName,
2274                                          theTimeStampId,
2275                                          &aNumDt,  
2276                                          &aNumOrd,
2277                                          &aDt);
2278         char profilename[MED_NAME_SIZE+1];
2279         med_int profilsize;
2280         char locname[MED_NAME_SIZE+1];
2281         MEDfieldnValueWithProfile(myFile->Id(),
2282                                   &aFieldName,
2283                                   aNumDt,
2284                                   aNumOrd,
2285                                   anEntity,
2286                                   med_geometry_type(aGeom),
2287                                   1,
2288                                   MED_COMPACT_STMODE,
2289                                   profilename,
2290                                   &profilsize,
2291                                   locname,
2292                                   &aNbGauss);
2293
2294         static TInt MAX_NB_GAUSS_POINTS = 32;
2295         if(aNbGauss <= 0 || aNbGauss > MAX_NB_GAUSS_POINTS)
2296           aNbGauss = 1;
2297
2298         aGeom2NbGauss[aGeom] = aNbGauss;
2299
2300         if(theErr) 
2301           *theErr = aRet;
2302         else if(aRet < 0)
2303           EXCEPTION(std::runtime_error,"GetTimeStampInfo - MEDfieldnValueWithProfile(...)");
2304       }      
2305     }
2306     
2307
2308     //----------------------------------------------------------------------------
2309     void 
2310     TVWrapper
2311     ::GetTimeStampValue(const PTimeStampValueBase& theTimeStampValue,
2312                         const TMKey2Profile& theMKey2Profile,
2313                         const TKey2Gauss& theKey2Gauss,
2314                         TErr* theErr)
2315     {
2316       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
2317       
2318       if(theErr && *theErr < 0)
2319         return;
2320       
2321       TIdt anId = myFile->Id();
2322       
2323       TValueHolder<EModeSwitch, med_switch_mode> aModeSwitch(theTimeStampValue->myModeSwitch);
2324       MED::TGeom2Profile& aGeom2Profile = theTimeStampValue->myGeom2Profile;
2325
2326       MED::PTimeStampInfo aTimeStampInfo = theTimeStampValue->myTimeStampInfo;
2327       TValueHolder<EEntiteMaillage, med_entity_type> anEntity(aTimeStampInfo->myEntity);
2328       TValueHolder<TInt, med_int> aNumDt(aTimeStampInfo->myNumDt);
2329       TValueHolder<TInt, med_int> aNumOrd(aTimeStampInfo->myNumOrd);
2330
2331       MED::PFieldInfo aFieldInfo = aTimeStampInfo->myFieldInfo;
2332       TValueHolder<TString, char> aFieldName(aFieldInfo->myName);
2333       TValueHolder<EBooleen, med_bool> anIsLocal(aFieldInfo->myIsLocal);
2334
2335       MED::PMeshInfo aMeshInfo = aFieldInfo->myMeshInfo;
2336       TValueHolder<TString, char> aMeshName(aMeshInfo->myName);
2337       
2338       TGeom2Gauss& aGeom2Gauss = aTimeStampInfo->myGeom2Gauss;
2339       TVector<char> aGaussName(GetNOMLength<eV2_2>()+1);
2340
2341       med_storage_mode aProfileMode = med_storage_mode(boost::get<0>(theMKey2Profile));
2342       MED::TKey2Profile aKey2Profile = boost::get<1>(theMKey2Profile);
2343       TVector<char> aProfileName(GetNOMLength<eV2_2>()+1);
2344
2345       TGeom2Size& aGeom2Size = aTimeStampInfo->myGeom2Size;
2346       TGeom2Size::iterator anIter = aGeom2Size.begin();
2347       for(; anIter != aGeom2Size.end(); anIter++){
2348         EGeometrieElement aGeom = anIter->first;
2349         TInt aNbElem = anIter->second;
2350         med_int profilesize,aNbGauss;
2351
2352         TInt aNbVal = MEDfieldnValueWithProfile(anId,
2353                                                 &aFieldName,
2354                                                 aNumDt,
2355                                                 aNumOrd,
2356                                                 anEntity,
2357                                                 med_geometry_type(aGeom),
2358                                                 1,
2359                                                 aProfileMode,
2360                                                 &aProfileName[0],
2361                                                 &profilesize,
2362                                                 &aGaussName[0],
2363                                                 &aNbGauss);
2364
2365         if(aNbVal <= 0){
2366           if(theErr){
2367             *theErr = -1;
2368             return;
2369           }
2370           EXCEPTION(std::runtime_error,"GetTimeStampValue - MEDfieldnValueWithProfile(...) - aNbVal == "<<aNbVal<<" <= 0");
2371         }
2372         
2373         TInt aNbComp = aFieldInfo->myNbComp;
2374         TInt aNbValue = aNbVal;// / aNbGauss; rules in MED changed
2375         theTimeStampValue->AllocateValue(aGeom,
2376                                          aNbValue,
2377                                          aNbGauss,
2378                                          aNbComp);
2379         TInt aValueSize = theTimeStampValue->GetValueSize(aGeom);
2380
2381         INITMSG(MYDEBUG,
2382                 "TVWrapper::GetTimeStampValue - aGeom = "<<aGeom<<
2383                 "; aNbVal = "<<aNbVal<<
2384                 "; aNbValue = "<<aNbValue<<
2385                 "; aNbGauss = "<<aNbGauss<<
2386                 "; aNbComp = "<<aNbComp<<
2387                 std::endl);
2388         
2389         TErr aRet = MEDfieldValueWithProfileRd(anId,
2390                                                &aFieldName,
2391                                                aNumDt,
2392                                                aNumOrd,
2393                                                anEntity,
2394                                                med_geometry_type(aGeom),
2395                                                aProfileMode,
2396                                                &aProfileName[0],
2397                                                aModeSwitch,
2398                                                MED_ALL_CONSTITUENT,
2399                                                theTimeStampValue->GetValuePtr(aGeom));
2400         if(aRet < 0){
2401           if(theErr){
2402             *theErr = MED_FALSE;
2403             return;
2404           }
2405           EXCEPTION(std::runtime_error,"GetTimeStampValue - MEDfieldValueWithProfileRd(...)");
2406         }
2407
2408         MED::PGaussInfo aGaussInfo;
2409         TGaussInfo::TKey aKey(aGeom,&aGaussName[0]);
2410         if(strcmp(&aGaussName[0],"") != 0){
2411           MED::TKey2Gauss::const_iterator anIter = theKey2Gauss.find(aKey);
2412           if(anIter != theKey2Gauss.end()){
2413             aGaussInfo = anIter->second;
2414             aGeom2Gauss[aGeom] = aGaussInfo;
2415           }
2416         }
2417         
2418         MED::PProfileInfo aProfileInfo;
2419         if(strcmp(&aProfileName[0],MED_NO_PROFILE) != 0){
2420           MED::TKey2Profile::const_iterator anIter = aKey2Profile.find(&aProfileName[0]);
2421           if(anIter != aKey2Profile.end()){
2422             aProfileInfo = anIter->second;
2423             aGeom2Profile[aGeom] = aProfileInfo;
2424           }
2425         }
2426
2427         if(aGaussInfo && aNbGauss != aGaussInfo->GetNbGauss()){
2428           if(theErr){
2429             *theErr = MED_FALSE;
2430             return;
2431           }
2432           EXCEPTION(std::runtime_error,"GetTimeStampValue - aNbGauss != aGaussInfo->GetNbGauss()");
2433         }
2434         
2435         if(aProfileInfo && aProfileInfo->IsPresent()){
2436           TInt aNbSubElem = aProfileInfo->GetSize();
2437           TInt aProfileSize = aNbSubElem*aNbComp*aNbGauss;
2438           if(aProfileSize != aValueSize){
2439             if(theErr){
2440               *theErr = -1;
2441               return;
2442             }
2443             EXCEPTION(std::runtime_error,
2444                       "GetTimeStampValue - aProfileSize("<<aProfileSize<<
2445                       ") != aValueSize("<<aValueSize<<
2446                       "); aNbVal = "<<aNbVal<<
2447                       "; anEntity = "<<anEntity<<
2448                       "; aGeom = "<<aGeom<<
2449                       "; aNbElem = "<<aNbElem<<
2450                       "; aNbSubElem = "<<aNbSubElem<<
2451                       "; aNbComp = "<<aNbComp<<
2452                       "; aNbGauss = "<<aNbGauss<<
2453                       "");
2454           }
2455         }else{
2456           if((aProfileMode == MED_GLOBAL_STMODE) && (aNbElem != aNbValue)){
2457             if(theErr){
2458               *theErr = -1;
2459               return;
2460             }
2461             EXCEPTION(std::runtime_error,
2462                       "GetTimeStampValue - aNbElem("<<aNbElem<<
2463                       ") != aNbValue("<<aNbValue<<
2464                       "); aNbVal = "<<aNbVal<<
2465                       "; anEntity = "<<anEntity<<
2466                       "; aGeom = "<<aGeom<<
2467                       "; aNbElem = "<<aNbElem<<
2468                       "; aNbComp = "<<aNbComp<<
2469                       "; aNbGauss = "<<aNbGauss<<
2470                       "");
2471           }
2472         }
2473       }
2474     }
2475     
2476     
2477     //----------------------------------------------------------------------------
2478     void
2479     TVWrapper
2480     ::SetTimeStampValue(const MED::PTimeStampValueBase& theTimeStampValue,
2481                         EModeAcces theMode,
2482                         TErr* theErr)
2483     {
2484       TFileWrapper aFileWrapper(myFile,theMode,theErr);
2485       
2486       if(theErr && *theErr < 0)
2487         return;
2488       
2489       TErr aRet;
2490       TIdt anId = myFile->Id();
2491       
2492       TValueHolder<EModeSwitch, med_switch_mode> aModeSwitch(theTimeStampValue->myModeSwitch);
2493       MED::TGeom2Profile& aGeom2Profile = theTimeStampValue->myGeom2Profile;
2494
2495       MED::PTimeStampInfo aTimeStampInfo = theTimeStampValue->myTimeStampInfo;
2496       TValueHolder<EEntiteMaillage, med_entity_type> anEntity(aTimeStampInfo->myEntity);
2497       TValueHolder<TInt, med_int> aNumDt(aTimeStampInfo->myNumDt);
2498       TValueHolder<TInt, med_int> aNumOrd(aTimeStampInfo->myNumOrd);
2499       TValueHolder<TString, char> anUnitDt(aTimeStampInfo->myUnitDt);
2500       TValueHolder<TFloat, med_float> aDt(aTimeStampInfo->myDt);
2501       MED::TGeom2Gauss& aGeom2Gauss = aTimeStampInfo->myGeom2Gauss;
2502
2503       MED::PFieldInfo aFieldInfo = aTimeStampInfo->myFieldInfo;
2504       TValueHolder<TString, char> aFieldName(aFieldInfo->myName);
2505
2506       MED::PMeshInfo aMeshInfo = aFieldInfo->myMeshInfo;
2507       TValueHolder<TString, char> aMeshName(aMeshInfo->myName);
2508       
2509       const TGeomSet& aGeomSet = theTimeStampValue->myGeomSet;
2510       TGeomSet::const_iterator anIter = aGeomSet.begin();
2511       for(; anIter != aGeomSet.end(); anIter++){
2512         EGeometrieElement aGeom = *anIter;
2513
2514         TVector<char> aGaussName(GetNOMLength<eV2_2>()+1);
2515         MED::TGeom2Gauss::const_iterator aGaussIter = aGeom2Gauss.find(aGeom);
2516         if(aGaussIter != aGeom2Gauss.end()){
2517           MED::PGaussInfo aGaussInfo = aGaussIter->second;
2518           strcpy(&aGaussName[0],&aGaussInfo->myName[0]);
2519         }
2520
2521         TVector<char> aProfileName(GetNOMLength<eV2_2>()+1);
2522         med_storage_mode aProfileMode = med_storage_mode(eNO_PFLMOD);
2523         MED::TGeom2Profile::const_iterator aProfileIter = aGeom2Profile.find(aGeom);
2524         if(aProfileIter != aGeom2Profile.end()){
2525           MED::PProfileInfo aProfileInfo = aProfileIter->second;
2526           aProfileMode = med_storage_mode(aProfileInfo->myMode);
2527           strcpy(&aProfileName[0],&aProfileInfo->myName[0]);
2528         }
2529
2530         med_int aNbVal = theTimeStampValue->GetNbVal(aGeom);
2531
2532         aRet = MEDfieldValueWithProfileWr(anId,
2533                                           &aFieldName,
2534                                           aNumDt,
2535                                           aNumOrd,
2536                                           aDt,
2537                                           anEntity,
2538                                           med_geometry_type(aGeom),
2539                                           aProfileMode,
2540                                           &aProfileName[0],
2541                                           &aGaussName[0],
2542                                           aModeSwitch,
2543                                           MED_ALL_CONSTITUENT,
2544                                           aNbVal,
2545                                           theTimeStampValue->GetValuePtr(aGeom));
2546         if(aRet < 0){
2547           if(theErr){
2548             *theErr = MED_FALSE;
2549             break;
2550           }
2551           EXCEPTION(std::runtime_error,"SetTimeStampValue - MEDfieldValueWithProfileWr(...)");
2552         }
2553         
2554       }
2555       
2556       INITMSG(MYDEBUG,"TVWrapper::SetTimeStampValue - MED_MODE_ACCES = "<<theMode<<"; aRet = "<<aRet<<std::endl);
2557     }
2558
2559     
2560     //----------------------------------------------------------------------------
2561     void 
2562     TVWrapper
2563     ::SetTimeStampValue(const PTimeStampValueBase& theTimeStampValue,
2564                         TErr* theErr)
2565     {
2566       TErr aRet;
2567       SetTimeStampValue(theTimeStampValue,eLECTURE_ECRITURE,&aRet);
2568       
2569       if(aRet < 0)
2570         SetTimeStampValue(theTimeStampValue,eLECTURE_AJOUT,&aRet);
2571
2572       if(theErr) 
2573         *theErr = aRet;
2574     }
2575     
2576     //----------------------------------------------------------------------------
2577     void 
2578     TVWrapper
2579     ::SetGrilleInfo(const MED::TGrilleInfo& theInfo,
2580                     TErr* theErr)
2581     {
2582       SetGrilleInfo(theInfo,eLECTURE_ECRITURE,theErr);
2583     }
2584
2585     //----------------------------------------------------------------------------
2586     void 
2587     TVWrapper
2588     ::SetGrilleInfo(const MED::TGrilleInfo& theInfo,
2589                     EModeAcces theMode,
2590                     TErr* theErr)
2591     {
2592       if(theInfo.myMeshInfo->myType != eSTRUCTURE)
2593         return;
2594       TFileWrapper aFileWrapper(myFile,theMode,theErr);
2595       
2596       if(theErr && *theErr < 0)
2597           return;
2598
2599       MED::TGrilleInfo& anInfo = const_cast<MED::TGrilleInfo&>(theInfo);
2600
2601       MED::TMeshInfo& aMeshInfo = *anInfo.myMeshInfo;
2602       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
2603
2604       TValueHolder<EGrilleType, med_grid_type > aGrilleType(anInfo.myGrilleType);
2605
2606       TErr aRet = 0;
2607       aRet = MEDmeshGridTypeRd(myFile->Id(),
2608                                &aMeshName,
2609                                &aGrilleType);
2610       if(theErr) 
2611         *theErr = aRet;
2612       else if(aRet < 0)
2613         EXCEPTION(std::runtime_error,"SetGrilleInfo - MEDmeshGridTypeRd(...)");
2614       
2615       if(anInfo.myGrilleType == eGRILLE_STANDARD){
2616         TValueHolder<TNodeCoord, med_float> aCoord(anInfo.myCoord);
2617         TValueHolder<EModeSwitch, med_switch_mode> aModeSwitch(anInfo.myModeSwitch);
2618         TValueHolder<TString, char> aCoordNames(anInfo.myCoordNames);
2619         TValueHolder<TString, char> aCoordUnits(anInfo.myCoordUnits);
2620         med_int aNbNoeuds = med_int(anInfo.myCoord.size() / aMeshInfo.myDim);
2621         //med_axis_type aRepere = MED_CARTESIAN;
2622
2623         aRet = MEDmeshNodeCoordinateWr(myFile->Id(),
2624                                        &aMeshName,
2625                                        MED_NO_DT,
2626                                        MED_NO_IT,
2627                                        MED_UNDEF_DT,
2628                                        aModeSwitch,
2629                                        aNbNoeuds,
2630                                        &aCoord);
2631
2632         if(aRet < 0)
2633           EXCEPTION(std::runtime_error,"SetGrilleInfo - MEDmeshNodeCoordinateWr(...)");
2634
2635         TValueHolder<TIntVector, med_int> aGrilleStructure(anInfo.myGrilleStructure);
2636         aRet = MEDmeshGridStructWr(myFile->Id(),
2637                                     &aMeshName,
2638                                    MED_NO_DT,
2639                                    MED_NO_IT,
2640                                    MED_UNDEF_DT,
2641                                    &aGrilleStructure);
2642         if(aRet < 0)
2643           EXCEPTION(std::runtime_error,"SetGrilleInfo - MEDmeshGridStructWr(...)");
2644         
2645       } else {
2646         for(med_int aAxis = 0; aAxis < aMeshInfo.myDim; aAxis++){
2647           aRet = MEDmeshGridIndexCoordinateWr(myFile->Id(),
2648                                               &aMeshName,
2649                                               MED_NO_DT,
2650                                               MED_NO_IT,
2651                                               MED_UNDEF_DT,
2652                                               aAxis+1,
2653                                               anInfo.GetIndexes(aAxis).size(),
2654                                               &anInfo.GetIndexes(aAxis)[0]);
2655
2656           if(aRet < 0)
2657             EXCEPTION(std::runtime_error,"SetGrilleInfo - MEDmeshGridIndexCoordinateWr(...)");
2658         }
2659         
2660       }
2661
2662       return;
2663     }
2664
2665     //----------------------------------------------------------------------------
2666     void
2667     TVWrapper
2668     ::GetGrilleInfo(TGrilleInfo& theInfo,
2669                     TErr* theErr)
2670     {
2671       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
2672
2673       if(theErr && *theErr < 0)
2674           return;
2675       
2676       MED::TMeshInfo& aMeshInfo = *theInfo.myMeshInfo;
2677       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
2678       EMaillage aMaillageType = aMeshInfo.myType;
2679       
2680       GetGrilleType(aMeshInfo, theInfo.myGrilleType, theErr);
2681       EGrilleType aGrilleType = theInfo.myGrilleType;
2682
2683       TErr aRet = 0;
2684       if(aMaillageType == eSTRUCTURE && aGrilleType == eGRILLE_STANDARD) {
2685         GetGrilleStruct(aMeshInfo, theInfo.myGrilleStructure, theErr);
2686
2687         TValueHolder<TNodeCoord, med_float> aCoord(theInfo.myCoord);
2688         TValueHolder<EModeSwitch, med_switch_mode> aModeSwitch(theInfo.myModeSwitch);
2689         TValueHolder<TString, char> aCoordNames(theInfo.myCoordNames);
2690         TValueHolder<TString, char> aCoordUnits(theInfo.myCoordUnits);
2691         //med_axis_type aRepere;
2692
2693         aRet = MEDmeshNodeCoordinateRd(myFile->Id(),
2694                                        &aMeshName,
2695                                        MED_NO_DT,
2696                                        MED_NO_IT,
2697                                        aModeSwitch,
2698                                        &aCoord);
2699
2700         if(theErr) 
2701           *theErr = aRet;
2702         else if(aRet < 0)
2703           EXCEPTION(std::runtime_error,"GetGrilleInfo - MEDmeshNodeCoordinateRd(...)");
2704
2705         //TInt aNbNodes = theInfo.GetNbNodes();//GetNbFamilies(aMeshInfo);
2706         TValueHolder<TElemNum, med_int> aFamNumNode(theInfo.myFamNumNode);
2707         
2708         aRet = MEDmeshEntityFamilyNumberRd(myFile->Id(),
2709                                            &aMeshName,
2710                                            MED_NO_DT,
2711                                            MED_NO_IT,
2712                                            MED_NODE,
2713                                            MED_NO_GEOTYPE,
2714                                            &aFamNumNode);
2715
2716         if(aRet < 0)
2717         {
2718 //            if (aRet == MED_ERR_DOESNTEXIST) // --- only valid with MED3.x files
2719               {
2720                 int mySize = (int)theInfo.myFamNumNode.size();
2721                 theInfo.myFamNumNode.clear();
2722                 theInfo.myFamNumNode.resize(mySize,0);
2723                 aRet = 0;
2724               }
2725 //            else
2726 //              EXCEPTION(std::runtime_error,"GetGrilleInfo - MEDmeshEntityFamilyNumberRd(...)");
2727         }
2728         if(theErr) 
2729           *theErr = aRet;
2730
2731         //============================
2732       }
2733
2734       if(aMaillageType == eSTRUCTURE && aGrilleType != eGRILLE_STANDARD){
2735         ETable aTable;
2736         for(med_int anAxis = 1; anAxis <= aMeshInfo.myDim; anAxis++){
2737           switch(anAxis){
2738           case 1 :
2739             aTable = eCOOR_IND1;
2740             break;
2741           case 2 :
2742             aTable = eCOOR_IND2;
2743             break;
2744           case 3 :
2745             aTable = eCOOR_IND3;
2746             break;
2747           default :
2748             aRet = -1;
2749           }
2750             
2751           if(theErr) 
2752             *theErr = aRet;
2753           else if(aRet < 0)
2754             EXCEPTION(std::runtime_error,"GetGrilleInfo - anAxis number out of range(...)");
2755           
2756           TInt aNbIndexes = GetNbNodes(aMeshInfo,aTable);
2757           if(aNbIndexes < 0)
2758             EXCEPTION(std::runtime_error,"GetGrilleInfo - Erreur a la lecture de la taille de l'indice");
2759             
2760           TValueHolder<TFloatVector, med_float> anIndexes(theInfo.GetIndexes(anAxis-1));
2761           //TValueHolder<ETable, med_data_type > table(aTable);
2762           //char aCompNames[MED_SNAME_SIZE+1];
2763           //char anUnitNames[MED_SNAME_SIZE+1];
2764           aRet=MEDmeshGridIndexCoordinateRd(myFile->Id(),&aMeshName,
2765                                             MED_NO_DT,MED_NO_IT,
2766                                             anAxis,
2767                                             &anIndexes);
2768
2769           //theInfo.SetCoordName(anAxis-1, aCompNames);
2770           //theInfo.SetCoordUnit(anAxis-1, anUnitNames);
2771           theInfo.SetGrilleStructure(anAxis-1, aNbIndexes);
2772
2773           if(theErr) 
2774             *theErr = aRet;
2775           else if(aRet < 0)
2776             EXCEPTION(std::runtime_error,"GetGrilleInfo - MEDindicesCoordLire(...)");
2777         }
2778       }
2779
2780       EGeometrieElement aGeom = theInfo.GetGeom();
2781       EEntiteMaillage aEntity = theInfo.GetEntity();
2782       TInt aNbCells = theInfo.GetNbCells();
2783       
2784       theInfo.myFamNum.resize(aNbCells);
2785       TValueHolder<TElemNum, med_int> aFamNum(theInfo.myFamNum);
2786       
2787       aRet = MEDmeshEntityFamilyNumberRd(myFile->Id(),
2788                                          &aMeshName,MED_NO_DT,MED_NO_IT,med_entity_type(aEntity),
2789                                          med_geometry_type(aGeom),&aFamNum);
2790
2791       if ( aMeshInfo.myDim == 3 )
2792       {
2793         aGeom = theInfo.GetSubGeom();
2794         aEntity = theInfo.GetSubEntity();
2795         aNbCells = theInfo.GetNbSubCells();
2796       
2797         theInfo.myFamSubNum.resize(aNbCells,0);
2798         TValueHolder<TElemNum, med_int> aFamNum(theInfo.myFamSubNum);
2799       
2800         aRet = MEDmeshEntityFamilyNumberRd(myFile->Id(),
2801                                     &aMeshName,MED_NO_DT,MED_NO_IT,
2802                                     med_entity_type(aEntity),
2803                                     med_geometry_type(aGeom),&aFamNum);
2804       }
2805       if(aRet < 0)
2806       {
2807 //          if (aRet == MED_ERR_DOESNTEXIST) // --- only valid with MED3.x files
2808             {
2809               int mySize = (int)theInfo.myFamNumNode.size();
2810               theInfo.myFamNumNode.clear();
2811               theInfo.myFamNumNode.resize(mySize,0);
2812               aRet = 0;
2813             }
2814 //          else
2815 //            EXCEPTION(std::runtime_error,"GetGrilleInfo - MEDmeshEntityFamilyNumberRd(...)");
2816       }
2817       if(theErr) 
2818         *theErr = aRet;
2819     }
2820
2821     void
2822     TVWrapper
2823     ::GetGrilleType(const MED::TMeshInfo& theMeshInfo,
2824                     EGrilleType& theGridType,
2825                     TErr* theErr)
2826     {
2827       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
2828
2829       if(theErr && *theErr < 0)
2830         EXCEPTION(std::runtime_error," GetGrilleType - aFileWrapper (...)");
2831
2832       MED::TMeshInfo& aMeshInfo = const_cast<MED::TMeshInfo&>(theMeshInfo);
2833       
2834       if(aMeshInfo.myType == eSTRUCTURE){
2835         TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
2836         TValueHolder<EGrilleType, med_grid_type> aGridType(theGridType);
2837         TErr aRet = MEDmeshGridTypeRd(myFile->Id(),
2838                                       &aMeshName,
2839                                       &aGridType);
2840
2841         if(aRet < 0)
2842           EXCEPTION(std::runtime_error,"GetGrilleInfo - MEDmeshGridTypeRd(...)");
2843       }
2844     }    
2845     
2846     void
2847     TVWrapper
2848     ::GetGrilleStruct(const MED::TMeshInfo& theMeshInfo,
2849                       TIntVector& theStruct,
2850                       TErr* theErr)
2851     {
2852       TFileWrapper aFileWrapper(myFile,eLECTURE,theErr);
2853       
2854       if(theErr && *theErr < 0)
2855           return;
2856       
2857       TErr aRet;
2858       MED::TMeshInfo& aMeshInfo = const_cast<MED::TMeshInfo&>(theMeshInfo);
2859
2860       TValueHolder<TString, char> aMeshName(aMeshInfo.myName);
2861       TValueHolder<TIntVector, med_int> aGridStructure(theStruct);
2862         
2863       aRet = MEDmeshGridStructRd(myFile->Id(),
2864                                  &aMeshName,
2865                                  MED_NO_DT,
2866                                  MED_NO_IT,
2867                                  &aGridStructure);
2868       if(theErr) 
2869         *theErr = aRet;
2870       else if(aRet < 0)
2871         EXCEPTION(std::runtime_error,"GetGrilleInfo - MEDmeshGridStructRd(...)");
2872     }
2873
2874   }  
2875 }