Salome HOME
Fix test case hang-up
[modules/paravis.git] / src / Plugins / MedReader / IO / vtkMedDriver30.cxx
1 // Copyright (C) 2010-2013  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "vtkMedDriver30.h"
21
22 #include "vtkObjectFactory.h"
23 #include "vtkStringArray.h"
24 #include "vtkDataArray.h"
25 #include "vtkIdTypeArray.h"
26 #include "vtkDoubleArray.h"
27 #include "vtkIntArray.h"
28 #include "vtkCharArray.h"
29
30 #include "vtkMedFile.h"
31 #include "vtkMedCartesianGrid.h"
32 #include "vtkMedPolarGrid.h"
33 #include "vtkMedCurvilinearGrid.h"
34 #include "vtkMedUnstructuredGrid.h"
35 #include "vtkMedField.h"
36 #include "vtkMedMesh.h"
37 #include "vtkMedFamily.h"
38 #include "vtkMedUtilities.h"
39 #include "vtkMedEntityArray.h"
40 #include "vtkMedLocalization.h"
41 #include "vtkMedProfile.h"
42 #include "vtkMedFieldOverEntity.h"
43 #include "vtkMedFieldStep.h"
44 #include "vtkMedGroup.h"
45 #include "vtkMedIntArray.h"
46 #include "vtkMedInterpolation.h"
47 #include "vtkMedFieldOnProfile.h"
48 #include "vtkMedFraction.h"
49 #include "vtkMedLink.h"
50 #include "vtkMedStructElement.h"
51 #include "vtkMedConstantAttribute.h"
52 #include "vtkMedVariableAttribute.h"
53
54 #include "vtkMultiProcessController.h"
55
56 #include <string>
57 #include <vector>
58 #include <algorithm>
59 using namespace std;
60
61 // vtkCxxRevisionMacro(vtkMedDriver30, "$Revision$")
62 vtkStandardNewMacro(vtkMedDriver30)
63
64 vtkMedDriver30::vtkMedDriver30() : vtkMedDriver()
65 {
66 }
67
68 vtkMedDriver30::~vtkMedDriver30()
69 {
70 }
71
72 // Description:
73 // load all Information data associated with this cartesian grid.
74 void vtkMedDriver30::ReadRegularGridInformation(vtkMedRegularGrid* grid)
75 {
76   FileOpen open(this);
77
78   grid->SetDimension(grid->GetParentMesh()->GetNumberOfAxis());
79
80   for(int axis=0; axis < grid->GetDimension(); axis++)
81     {
82     med_int size;
83     med_bool coordinatechangement, geotransformation;
84
85     if ((size = MEDmeshnEntity(
86                   this->FileId,
87                   grid->GetParentMesh()->GetName(),
88                   grid->GetComputeStep().TimeIt,
89                   grid->GetComputeStep().IterationIt,
90                   MED_NODE,
91                   MED_NONE,
92                   (med_data_type)(MED_COORDINATE_AXIS1 + axis),
93                   MED_NO_CMODE,
94                   &coordinatechangement,
95                   &geotransformation)) < 0)
96       {
97       vtkErrorMacro("ERROR : number of coordinates on X axis ...");
98       }
99
100     grid->SetAxisSize(axis, size);
101     }
102
103   med_int ncell = 1;
104   if(grid->GetAxisSize(0) > 1)
105     ncell *= grid->GetAxisSize(0)-1;
106   if(grid->GetAxisSize(1) > 1)
107     ncell *= grid->GetAxisSize(1)-1;
108   if(grid->GetAxisSize(2) > 1)
109     ncell *= grid->GetAxisSize(2)-1;
110
111   vtkMedEntity entity;
112   entity.EntityType = MED_CELL;
113
114   switch(grid->GetDimension())
115     {
116     case 0 :
117       entity.GeometryType = MED_POINT1;
118       break;
119     case 1 :
120       entity.GeometryType = MED_SEG2;
121       break;
122     case 2 :
123       entity.GeometryType = MED_QUAD4;
124       break;
125     case 3 :
126       entity.GeometryType = MED_HEXA8;
127       break;
128     default :
129         vtkErrorMacro("Unsupported dimension for curvilinear grid : "
130                       << grid->GetDimension());
131     return;
132     }
133
134   vtkMedEntityArray* array = vtkMedEntityArray::New();
135   array->SetParentGrid(grid);
136   array->SetNumberOfEntity(ncell);
137   array->SetEntity(entity);
138   array->SetConnectivity(MED_NODAL);
139   grid->AppendEntityArray(array);
140   array->Delete();
141   // this triggers the creation of undefined families
142   this->LoadFamilyIds(array);
143
144 }
145
146 void  vtkMedDriver30::LoadRegularGridCoordinates(vtkMedRegularGrid* grid)
147 {
148   FileOpen open(this);
149
150   for(int axis=0; axis < grid->GetParentMesh()->GetNumberOfAxis(); axis++)
151     {
152
153     vtkDataArray* coords = vtkMedUtilities::NewCoordArray();
154     grid->SetAxisCoordinate(axis, coords);
155     coords->Delete();
156
157     coords->SetNumberOfComponents(1);
158     coords->SetNumberOfTuples(grid->GetAxisSize(axis));
159
160     if (MEDmeshGridIndexCoordinateRd(
161           this->FileId,
162           grid->GetParentMesh()->GetName(),
163           grid->GetComputeStep().TimeIt,
164           grid->GetComputeStep().IterationIt,
165           axis+1,
166           (med_float*)coords->GetVoidPointer(0)) < 0)
167       {
168       vtkErrorMacro("ERROR : read axis " << axis << " coordinates ...");
169       grid->SetAxisCoordinate(axis, NULL);
170       return;
171       }
172     }
173
174 }
175
176 // Description:
177 // load all Information data associated with this standard grid.
178 void vtkMedDriver30::ReadCurvilinearGridInformation(vtkMedCurvilinearGrid* grid)
179 {
180   FileOpen open(this);
181
182   grid->SetDimension(grid->GetParentMesh()->GetNumberOfAxis());
183
184   med_int size;
185   med_bool coordinatechangement, geotransformation;
186
187   if ((size = MEDmeshnEntity(
188                 this->FileId,
189                 grid->GetParentMesh()->GetName(),
190                 grid->GetComputeStep().TimeIt,
191                 grid->GetComputeStep().IterationIt,
192                 MED_NODE,
193                 MED_NONE,
194                 MED_COORDINATE,
195                 MED_NO_CMODE,
196                 &coordinatechangement,
197                 &geotransformation)) < 0)
198     {
199     vtkErrorMacro("ReadCurvilinearGridInformation MEDmeshnEntity");
200     }
201
202   grid->SetNumberOfPoints(size);
203
204   med_int axissize[3];
205
206   if(MEDmeshGridStructRd(
207       this->FileId,
208       grid->GetParentMesh()->GetName(),
209       grid->GetComputeStep().TimeIt,
210       grid->GetComputeStep().IterationIt,
211       axissize) < 0)
212     {
213     vtkErrorMacro("ReadCurvilinearGridInformation MEDmeshGridStructRd");
214     }
215
216   switch(grid->GetDimension())
217     {
218     case 0 :
219       break;
220     case 1 :
221       grid->SetAxisSize(0, (axissize[0] <= 0 ? 1: axissize[0]));
222       break;
223     case 2 :
224       grid->SetAxisSize(0, (axissize[0] <= 0 ? 1: axissize[0]));
225       grid->SetAxisSize(1, (axissize[1] <= 0 ? 1: axissize[1]));
226       break;
227     case 3 :
228       grid->SetAxisSize(0, (axissize[0] <= 0 ? 1: axissize[0]));
229       grid->SetAxisSize(1, (axissize[1] <= 0 ? 1: axissize[1]));
230       grid->SetAxisSize(2, (axissize[2] <= 0 ? 1: axissize[2]));
231       break;
232     default :
233         vtkErrorMacro("Unsupported dimension for curvilinear grid : "
234                       << grid->GetDimension());
235     return;
236     }
237
238   // A test to verify the number of points : total number of
239   // points should be equal to the product of each axis size
240   med_int size2 = 0;
241
242   if (grid->GetDimension() == 1)
243     {
244     size2 = grid->GetAxisSize(0);
245     }
246
247   if (grid->GetDimension() == 2)
248     {
249     size2 = grid->GetAxisSize(0)*grid->GetAxisSize(1);
250     }
251
252   if (grid->GetDimension() == 3)
253     {
254     size2 = grid->GetAxisSize(0)*grid->GetAxisSize(1)*grid->GetAxisSize(2);
255     }
256
257   if(size != size2)
258     {
259     vtkErrorMacro("The total number of points of a Curvilinear grid should "
260                   << "be the product of each axis size!");
261     }
262
263   med_int ncell = 1;
264   if ((grid->GetDimension() >= 1) && (grid->GetAxisSize(0) > 1))
265     ncell *= grid->GetAxisSize(0)-1;
266   if ((grid->GetDimension() >= 2) && (grid->GetAxisSize(1) > 1))
267     ncell *= grid->GetAxisSize(1)-1;
268   if ((grid->GetDimension() >= 3) && (grid->GetAxisSize(2) > 1))
269     ncell *= grid->GetAxisSize(2)-1;
270
271   vtkMedEntity entity;
272   entity.EntityType = MED_CELL;
273
274   switch(grid->GetDimension())
275     {
276     case 0 :
277       entity.GeometryType = MED_POINT1;
278       break;
279     case 1 :
280       entity.GeometryType = MED_SEG2;
281       break;
282     case 2 :
283       entity.GeometryType = MED_QUAD4;
284       break;
285     case 3 :
286       entity.GeometryType = MED_HEXA8;
287       break;
288     default :
289         vtkErrorMacro("Unsupported dimension for curvilinear grid : "
290                       << grid->GetDimension());
291     return;
292     }
293
294   vtkMedEntityArray* array = vtkMedEntityArray::New();
295   array->SetParentGrid(grid);
296   array->SetNumberOfEntity(ncell);
297   array->SetEntity(entity);
298   array->SetConnectivity(MED_NODAL);
299   grid->AppendEntityArray(array);
300   array->Delete();
301   // this triggers the creation of undefined families
302   this->LoadFamilyIds(array);
303 }
304
305 // Description : read the number of entity of all geometry type
306 // for a given entity type and a given connectivity mode
307 void vtkMedDriver30::ReadNumberOfEntity(
308     vtkMedUnstructuredGrid* grid,
309     med_entity_type entityType,
310     med_connectivity_mode connectivity)
311 {
312   FileOpen open(this);
313
314   med_bool changement, transformation;
315
316   const char* meshName = grid->GetParentMesh()->GetName();
317
318   const vtkMedComputeStep& cs = grid->GetComputeStep();
319
320   med_int nentity = MEDmeshnEntity(
321                         this->FileId,
322                         meshName,
323                         cs.TimeIt,
324                         cs.IterationIt,
325                         entityType,
326                         MED_GEO_ALL,
327                         MED_UNDEF_DATATYPE ,
328                         connectivity,
329                         &changement,
330                         &transformation );
331
332   for(med_int geotypeit = 1; geotypeit <= nentity; geotypeit++)
333     {
334     // read cell informations
335     vtkMedEntity entity;
336     entity.EntityType = entityType;
337
338     char geometryName[MED_NAME_SIZE+1] = "";
339
340     // this gives us the med_geometry_type
341     if( MEDmeshEntityInfo( FileId, meshName,
342                            cs.TimeIt,
343                            cs.IterationIt,
344                            entityType,
345                            geotypeit,
346                            geometryName,
347                            &entity.GeometryType) < 0)
348       {
349       vtkErrorMacro("MEDmeshEntityInfo");
350       continue;
351       }
352
353     entity.GeometryName = geometryName;
354     med_int ncell = 0;
355
356     if(entity.GeometryType == MED_POLYGON)
357       {
358       // read the number of cells of this type
359       ncell = MEDmeshnEntity(this->FileId,
360                              meshName,
361                              cs.TimeIt,
362                              cs.IterationIt,
363                              entity.EntityType,
364                              entity.GeometryType,
365                              MED_INDEX_NODE,
366                              connectivity,
367                              &changement,
368                              &transformation ) - 1;
369       }
370     else if(entity.GeometryType == MED_POLYHEDRON)
371       {
372       // read the number of cells of this type
373       ncell = MEDmeshnEntity(this->FileId,
374                              meshName,
375                              cs.TimeIt,
376                              cs.IterationIt,
377                              entity.EntityType,
378                              entity.GeometryType,
379                              MED_INDEX_FACE,
380                              connectivity,
381                              &changement,
382                              &transformation  ) - 1;
383       }
384     else
385       {
386       ncell = MEDmeshnEntity(this->FileId,
387                              meshName,
388                              cs.TimeIt,
389                              cs.IterationIt,
390                              entity.EntityType,
391                              entity.GeometryType,
392                              MED_CONNECTIVITY,
393                              connectivity,
394                              &changement,
395                              &transformation  );
396       }
397
398     if(ncell > 0)
399       {
400       vtkMedEntityArray* array = vtkMedEntityArray::New();
401       array->SetParentGrid(grid);
402       array->SetNumberOfEntity(ncell);
403       array->SetEntity(entity);
404       array->SetConnectivity(connectivity);
405       grid->AppendEntityArray(array);
406       array->Delete();
407       // this triggers the creation of undefined families
408       this->LoadFamilyIds(array);
409       }
410     }
411 }
412
413 // Description:
414 // load all Information data associated with this unstructured grid.
415 void vtkMedDriver30::ReadUnstructuredGridInformation(
416     vtkMedUnstructuredGrid* grid)
417 {
418   FileOpen open(this);
419
420   vtkMedMesh *mesh = grid->GetParentMesh();
421
422   const char *meshName = mesh->GetName();
423   med_connectivity_mode connectivity;
424
425   med_bool changement;
426   med_bool transformation;
427   med_int profilesize;
428
429   char profilename[MED_NAME_SIZE+1];
430   memset(profilename, '\0', MED_NAME_SIZE+1);
431
432   vtkMedComputeStep cs = grid->GetComputeStep();
433
434   // first check if we have points
435   vtkIdType npoints = MEDmeshnEntityWithProfile(
436                         this->FileId,
437                         meshName,
438                         cs.TimeIt,
439                         cs.IterationIt,
440                         MED_NODE,
441                         MED_NONE,
442                         MED_COORDINATE,
443                         MED_NODAL,
444                         MED_COMPACT_PFLMODE,
445                         profilename,
446                         &profilesize,
447                         &changement,
448                         &transformation);
449
450   if(npoints > 0)
451     {
452     grid->SetNumberOfPoints(npoints);
453     }
454   else
455     {
456     if(grid->GetPreviousGrid() == NULL)
457       {
458       vtkErrorMacro("No point and no previous grid");
459       }
460     grid->SetUsePreviousCoordinates(true);
461     grid->SetNumberOfPoints(grid->GetPreviousGrid()->GetNumberOfPoints());
462     return;
463     }
464
465   this->ReadNumberOfEntity(grid, MED_CELL, MED_NODAL);
466   //this->ReadNumberOfEntity(grid, MED_CELL, MED_DESCENDING);  // VSR: 29.11.2013: commented according to the issue 22327: descending connectivity is wrongly processed
467   this->ReadNumberOfEntity(grid, MED_DESCENDING_FACE, MED_NODAL);
468   this->ReadNumberOfEntity(grid, MED_DESCENDING_FACE, MED_DESCENDING);
469   this->ReadNumberOfEntity(grid, MED_DESCENDING_EDGE, MED_NODAL);
470   this->ReadNumberOfEntity(grid, MED_DESCENDING_EDGE, MED_DESCENDING);
471   this->ReadNumberOfEntity(grid, MED_STRUCT_ELEMENT, MED_NODAL);
472   this->ReadNumberOfEntity(grid, MED_STRUCT_ELEMENT, MED_DESCENDING);
473
474   // create the point vtkMedEntityArray support
475   vtkMedEntity entity;
476   entity.EntityType = MED_NODE;
477   entity.GeometryType = MED_POINT1;
478   vtkMedEntityArray* pea = vtkMedEntityArray::New();
479   pea->SetEntity(entity);
480   pea->SetParentGrid(grid);
481   pea->SetNumberOfEntity(grid->GetNumberOfPoints());
482   grid->AppendEntityArray(pea);
483   pea->Delete();
484
485   this->LoadFamilyIds(pea);
486 }
487
488 void vtkMedDriver30::ReadFamilyInformation(vtkMedMesh* mesh, vtkMedFamily* family)
489 {
490   FileOpen open(this);
491
492   med_int familyid;
493   med_int ngroup;
494   char* groupNames = NULL;
495   const  char* meshName = mesh->GetName();
496
497   ngroup = MEDnFamilyGroup(FileId, meshName, family->GetMedIterator());
498
499   bool has_no_group = false;
500   if(ngroup <= 0)
501     {
502     if(ngroup < 0)
503       {
504       vtkErrorMacro("Error while reading the number of groups");
505       }
506     ngroup = 1;
507     has_no_group = true;
508     }
509
510   groupNames = new char[ngroup * MED_LNAME_SIZE + 1];
511   memset(groupNames, '\0', ngroup * MED_LNAME_SIZE + 1);
512
513   // special case for files written by med < 3,
514   // I have to use the 23 interface
515   if(mesh->GetParentFile()->GetVersionMajor() < 3)
516     {
517     med_int *attributenumber = NULL;
518     med_int *attributevalue = NULL;
519     char *attributedes = NULL;
520
521     med_int nattr = MEDnFamily23Attribute(
522                       this->FileId,
523                       meshName,
524                       family->GetMedIterator());
525
526     if(nattr < 0)
527       {
528       vtkErrorMacro("MEDnFamily23Attribute");
529       }
530
531     if(nattr > 0)
532       {
533       attributenumber = new med_int[nattr];
534       attributevalue = new med_int[nattr];
535       attributedes = new char[nattr*MED_COMMENT_SIZE+1];
536       memset(attributedes, '\0', nattr*MED_COMMENT_SIZE+1);
537       }
538
539     char familyName[MED_NAME_SIZE+1] = "";
540
541     if(MEDfamily23Info (this->FileId,
542                         meshName,
543                         family->GetMedIterator(),
544                         familyName,
545                         attributenumber,
546                         attributevalue,
547                         attributedes,
548                         &familyid,
549                         groupNames ) < 0)
550       {
551       vtkDebugMacro("MEDfamily23Info");
552       }
553
554     family->SetName(familyName);
555
556     if(attributenumber != NULL)
557       delete[] attributenumber;
558     if(attributevalue != NULL)
559       delete[] attributevalue;
560     if(attributedes != NULL)
561       delete[] attributedes;
562     }
563   else
564     {
565     char familyName[MED_NAME_SIZE+1] = "";
566     if(MEDfamilyInfo( this->FileId,
567                       meshName,
568                       family->GetMedIterator(),
569                       familyName,
570                       &familyid,
571                       groupNames ) < 0)
572       {
573       vtkErrorMacro(
574           "vtkMedDriver30::ReadInformation(vtkMedFamily* family)"
575           << " cannot read family informations.");
576       return;
577       }
578     family->SetName(familyName);
579     }
580
581   family->SetId(familyid);
582
583   // Fix for the issue "0021721: [CEA 590] Invalid groups on cells
584   if( familyid == 0 ) {
585     family->SetPointOrCell( vtkMedUtilities::OnCell );
586     mesh->AppendCellFamily( family );
587   } else {
588     //rnv: improve algorithm to determine entity of the family:
589     //     1) Read Nb nodes
590     //     2) Read families of the nodes
591     //     3) If list of the families of the nodes contains familyid => vtkMedUtilities::OnPoint
592     //     otherwise => vtkMedUtilities::OnCell
593     med_bool  v1, v2;
594     med_int size = MEDmeshnEntity(
595                                   this->FileId,
596                                   meshName,
597                                   MED_NO_DT,
598                                   MED_NO_IT,
599                                   MED_NODE,
600                                   MED_NO_GEOTYPE,
601                                   MED_COORDINATE,
602                                   MED_NO_CMODE,
603                                   &v1,
604                                   &v2 );
605     if( size < 0 ) {
606       vtkErrorMacro( "vtkMedDriver30::ReadInformation(vtkMedFamily* family)"
607                      <<" cannot read nb Nodes" );
608       return;
609     }
610     
611     vector<med_int> n_fams;
612     n_fams.resize( size );
613     
614     med_int ret_val = MEDmeshEntityFamilyNumberRd( this->FileId,
615                                                    meshName,
616                                                    MED_NO_DT,
617                                                    MED_NO_IT,
618                                                    MED_NODE,
619                                                    MED_NO_GEOTYPE ,
620                                                    &n_fams[0] );
621     // Remove ZERO FAMILY
622     remove( n_fams.begin(),n_fams.end(), 0 );
623
624     bool isOnPoints = ( ret_val >= 0) && (find( n_fams.begin(), n_fams.end(), familyid ) != n_fams.end() );
625     if( isOnPoints ) {
626       family->SetPointOrCell(vtkMedUtilities::OnPoint);
627       mesh->AppendPointFamily(family);
628     } else {
629       family->SetPointOrCell(vtkMedUtilities::OnCell);        
630       mesh->AppendCellFamily(family);   
631     }
632   }
633   
634   family->AllocateNumberOfGroup(ngroup);
635   // if there where no group, set the name to the default value
636   if(has_no_group)
637     {
638     memcpy(groupNames, vtkMedUtilities::NoGroupName,
639            strlen(vtkMedUtilities::NoGroupName));
640     }
641
642   for(int index = 0; index < ngroup; index++)
643     {
644     char realGroupName[MED_LNAME_SIZE + 1];
645     memset(realGroupName, '\0', MED_LNAME_SIZE + 1);
646     memcpy(realGroupName, groupNames + index * MED_LNAME_SIZE,
647         MED_LNAME_SIZE * sizeof(char));
648     vtkMedGroup* group = mesh->GetOrCreateGroup(family->GetPointOrCell(),
649         realGroupName);
650
651     family->SetGroup(index, group);
652     }
653
654   delete[] groupNames;
655
656   if(familyid == 0)
657     {
658     vtkMedFamily* famzero = vtkMedFamily::New();
659     mesh->AppendPointFamily(famzero);
660     famzero->Delete();
661
662     famzero->SetName(family->GetName());
663     famzero->SetMedIterator(family->GetMedIterator());
664     famzero->SetId(family->GetId());
665     famzero->SetPointOrCell(vtkMedUtilities::OnPoint);
666     famzero->AllocateNumberOfGroup(family->GetNumberOfGroup());
667     for(int gid=0; gid<family->GetNumberOfGroup(); gid++)
668       {
669       vtkMedGroup* group = mesh->GetOrCreateGroup(
670           vtkMedUtilities::OnPoint,
671           family->GetGroup(gid)->GetName());
672       famzero->SetGroup(gid, group);
673       mesh->AppendPointGroup(group);
674       }
675     }
676 }
677
678 void  vtkMedDriver30::ReadLinkInformation(vtkMedLink* link)
679 {
680   med_int size;
681   char linkMeshName[MED_NAME_SIZE+1] = "";
682   if(MEDlinkInfo(this->FileId,
683                  link->GetMedIterator(),
684                  linkMeshName,
685                  &size) < 0)
686     {
687     vtkErrorMacro("MEDlinkInfo");
688     return;
689     }
690   link->SetMeshName(linkMeshName);
691   if(size <= 0)
692     return;
693
694   char* path = new char[size + 1];
695   memset(path, '\0', size+1);
696   if(MEDlinkRd(this->FileId, link->GetMeshName(), path) < 0)
697     {
698     vtkErrorMacro("MEDlinkRd");
699     memset(path, '\0', size+1);
700     }
701
702   link->SetLink(path);
703
704   delete[] path;
705 }
706
707 void vtkMedDriver30::ReadFileInformation(vtkMedFile* file)
708 {
709   FileOpen open(this);
710
711   char comment[MED_COMMENT_SIZE+1] = "";
712
713   MEDfileCommentRd(this->FileId,
714                   comment);
715
716   file->SetComment(comment);
717
718   med_int major, minor, release;
719   MEDfileNumVersionRd(this->FileId, &major, &minor, &release);
720   file->SetVersionMajor(major);
721   file->SetVersionMinor(minor);
722   file->SetVersionRelease(release);
723
724   int nlink = MEDnLink(this->FileId);
725   file->AllocateNumberOfLink(nlink);
726   for(int linkid=0; linkid<nlink; linkid++)
727     {
728     vtkMedLink* link = file->GetLink(linkid);
729     link->SetMedIterator(linkid+1);
730     this->ReadLinkInformation(link);
731     }
732
733   int nprof = MEDnProfile(FileId);
734   // Reading id s not possible in parallel if the file contains Profiles
735   vtkMultiProcessController* controller = vtkMultiProcessController::GetGlobalController();
736   if (controller != NULL)
737     if ((nprof != 0) && (controller->GetNumberOfProcesses() > 1))
738     {
739       vtkWarningMacro("ATTENTION: The MED Reader cannot read profiles when used in parallel");
740     return;
741     }
742   file->AllocateNumberOfProfile(nprof);
743   for(int i = 0; i < nprof; i++)
744     {
745     vtkMedProfile* profile = file->GetProfile(i);
746     profile->SetMedIterator(i + 1);
747     profile->SetParentFile(file);
748     this->ReadProfileInformation(profile);
749     }
750
751   int nloc = MEDnLocalization(this->FileId);
752   file->AllocateNumberOfLocalization(nloc);
753   for(int i = 0; i < nloc; i++)
754     {
755     vtkMedLocalization* loc = file->GetLocalization(i);
756     loc->SetMedIterator(i + 1);
757     loc->SetParentFile(file);
758     this->ReadLocalizationInformation(loc);
759     }
760
761   int nsupportmesh = MEDnSupportMesh(this->FileId);
762   file->AllocateNumberOfSupportMesh(nsupportmesh);
763   for(int i = 0; i < nsupportmesh; i++)
764     {
765     vtkMedMesh* supportmesh = file->GetSupportMesh(i);
766     supportmesh->SetMedIterator(i + 1);
767     supportmesh->SetParentFile(file);
768     this->ReadSupportMeshInformation(supportmesh);
769     }
770
771   int nmesh = MEDnMesh(this->FileId);
772   file->AllocateNumberOfMesh(nmesh);
773   for(int i = 0; i < nmesh; i++)
774     {
775     vtkMedMesh* mesh = file->GetMesh(i);
776     mesh->SetMedIterator(i + 1);
777     mesh->SetParentFile(file);
778     this->ReadMeshInformation(mesh);
779     }
780
781   int nfields = MEDnField(this->FileId);
782   file->AllocateNumberOfField(nfields);
783   for(int i = 0; i < nfields; i++)
784     {
785     vtkMedField* field = file->GetField(i);
786     field->SetMedIterator(i + 1);
787     field->SetParentFile(file);
788     this->ReadFieldInformation(field);
789     field->ComputeFieldType();
790     while(field->HasManyFieldTypes())
791       {
792       vtkMedField* newfield = vtkMedField::New();
793       int type = field->GetFirstType();
794       newfield->ExtractFieldType(field, type);
795       file->AppendField(newfield);
796       newfield->Delete();
797       }
798     }
799
800   int nstruct = MEDnStructElement(this->FileId);
801
802   file->AllocateNumberOfStructElement(nstruct);
803   for(int i = 0; i < nstruct; i++)
804     {
805     vtkMedStructElement* structelem = file->GetStructElement(i);
806     structelem->SetMedIterator(i+1);
807     structelem->SetParentFile(file);
808     this->ReadStructElementInformation(structelem);
809     }
810
811 }
812
813 void  vtkMedDriver30::ReadStructElementInformation(
814     vtkMedStructElement* structelem)
815 {
816
817   FileOpen open(this);
818
819   char modelname[MED_NAME_SIZE+1] = "";
820   med_geometry_type mgeotype;
821   med_int modeldim;
822   char supportmeshname[MED_NAME_SIZE+1] = "";
823   med_entity_type sentitytype;
824   med_int snbofnode;
825   med_int snbofcell;
826   med_geometry_type sgeotype;
827   med_int nbofconstantattribute;
828   med_bool anyprofile;
829   med_int nbofvariableattribute;
830
831   if(MEDstructElementInfo (this->FileId,
832                            structelem->GetMedIterator(),
833                            modelname,
834                            &mgeotype,
835                            &modeldim,
836                            supportmeshname,
837                            &sentitytype,
838                            &snbofnode,
839                            &snbofcell,
840                            &sgeotype,
841                            &nbofconstantattribute,
842                            &anyprofile,
843                            &nbofvariableattribute) < 0)
844     {
845     vtkErrorMacro("Error in MEDstructElementInfo");
846     return;
847     }
848   structelem->SetName(modelname);
849   structelem->SetGeometryType(mgeotype);
850   structelem->SetModelDimension(modeldim);
851   structelem->SetSupportMeshName(supportmeshname);
852   structelem->SetSupportEntityType(sentitytype);
853   structelem->SetSupportNumberOfNode(snbofnode);
854   structelem->SetSupportNumberOfCell(snbofcell);
855   structelem->SetSupportGeometryType(sgeotype);
856   structelem->AllocateNumberOfConstantAttribute(nbofconstantattribute);
857   structelem->AllocateNumberOfVariableAttribute(nbofvariableattribute);
858   structelem->SetAnyProfile(anyprofile);
859
860   for(int attit = 0; attit < nbofconstantattribute; attit ++)
861     {
862     vtkMedConstantAttribute* constatt = structelem->GetConstantAttribute(attit);
863     constatt->SetMedIterator(attit+1);
864     constatt->SetParentStructElement(structelem);
865     this->ReadConstantAttributeInformation(constatt);
866     }
867
868   for(int attit = 0; attit < nbofvariableattribute; attit ++)
869     {
870     vtkMedVariableAttribute* varatt = structelem->GetVariableAttribute(attit);
871     varatt->SetMedIterator(attit+1);
872     varatt->SetParentStructElement(structelem);
873     this->ReadVariableAttributeInformation(varatt);
874     }
875 }
876
877 void vtkMedDriver30::ReadConstantAttributeInformation(
878     vtkMedConstantAttribute* constAttr)
879 {
880
881   FileOpen open(this);
882
883   char constattname[MED_NAME_SIZE+1] = "";
884   med_attribute_type constatttype;
885   med_int nbofcomponent;
886   med_entity_type sentitytype;
887   char profilename[MED_NAME_SIZE+1] = "";
888   med_int profilesize;
889
890   if(MEDstructElementConstAttInfo(
891       this->FileId,
892       constAttr->GetParentStructElement()->GetName(),
893       constAttr->GetMedIterator(),
894       constattname,
895       &constatttype,
896       &nbofcomponent,
897       &sentitytype,
898       profilename,
899       &profilesize)   < 0)
900     {
901     vtkErrorMacro("MEDstructElementConstAttInfo error");
902     return;
903     }
904
905   constAttr->SetName(constattname);
906   constAttr->SetAttributeType(constatttype);
907   constAttr->SetNumberOfComponent(nbofcomponent);
908   constAttr->SetSupportEntityType(sentitytype);
909   constAttr->SetProfileName(profilename);
910   constAttr->SetProfileSize(profilesize);
911
912   vtkAbstractArray* values = vtkMedUtilities::NewArray(constatttype);
913   if(values == NULL)
914     return;
915   constAttr->SetValues(values);
916   values->Delete();
917
918   values->SetNumberOfComponents(nbofcomponent);
919   vtkIdType ntuple = 0;
920   if((strcmp(profilename, MED_NO_PROFILE) != 0) &&
921      (strcmp(profilename, "\0") != 0))
922     {
923     ntuple = profilesize;
924     }
925   else if(constAttr->GetSupportEntityType() == MED_CELL)
926     {
927     ntuple = constAttr->GetParentStructElement()->GetSupportNumberOfCell();
928     }
929   else
930     {
931     ntuple = constAttr->GetParentStructElement()->GetSupportNumberOfNode();
932     }
933   values->SetNumberOfTuples(ntuple);
934
935   void* ptr = NULL;
936   vtkSmartPointer<vtkCharArray> buffer = vtkSmartPointer<vtkCharArray>::New();
937   if(constatttype != MED_ATT_NAME)
938     {
939     ptr = values->GetVoidPointer(0);
940     }
941   else
942     {
943     buffer->SetNumberOfValues(MED_NAME_SIZE*nbofcomponent*ntuple);
944     ptr = buffer->GetVoidPointer(0);
945     }
946
947   if(MEDstructElementConstAttRd (this->FileId,
948         constAttr->GetParentStructElement()->GetName(),
949         constAttr->GetName(), ptr) < 0)
950     {
951     vtkErrorMacro("MEDstructElementConstAttRd");
952     return;
953     }
954
955   if(constatttype == MED_ATT_NAME)
956     {
957     char name[MED_NAME_SIZE+1] = "";
958     char* nameptr = (char*) ptr;
959     vtkStringArray* names = vtkStringArray::SafeDownCast(values);
960     for(vtkIdType id = 0; id < nbofcomponent*ntuple; id++)
961       {
962       memset(name, '\0', MED_NAME_SIZE+1);
963       strncpy(name, nameptr + id * MED_NAME_SIZE, MED_NAME_SIZE);
964       names->SetValue(id, name);
965       }
966     }
967
968   return;
969 }
970
971 void vtkMedDriver30::ReadVariableAttributeInformation(
972     vtkMedVariableAttribute* varAttr)
973 {
974
975   FileOpen open(this);
976
977   char varattname[MED_NAME_SIZE+1] = "";
978   med_attribute_type varatttype;
979   med_int nbofcomponent;
980
981   if(MEDstructElementVarAttInfo (
982       this->FileId,
983       varAttr->GetParentStructElement()->GetName(),
984       varAttr->GetMedIterator(),
985       varattname,
986       &varatttype,
987       &nbofcomponent) < 0)
988     {
989     vtkErrorMacro("MEDstructElementVarAttInfo");
990     return;
991     }
992
993   varAttr->SetName(varattname);
994   varAttr->SetAttributeType(varatttype);
995   varAttr->SetNumberOfComponent(nbofcomponent);
996
997   return;
998 }
999
1000 void vtkMedDriver30::ReadProfileInformation(vtkMedProfile* profile)
1001 {
1002   FileOpen open(this);
1003
1004   med_int nelem;
1005   char profileName[MED_NAME_SIZE+1] = "";
1006
1007   if(MEDprofileInfo(this->FileId,
1008                     profile->GetMedIterator(),
1009                     profileName,
1010                     &nelem) < 0)
1011     {
1012     vtkErrorMacro("cannot read information on profile"
1013         << profile->GetMedIterator());
1014     }
1015   profile->SetName(profileName);
1016   profile->SetNumberOfElement(nelem);
1017 }
1018
1019 void vtkMedDriver30::ReadFieldInformation(vtkMedField* field)
1020 {
1021   FileOpen open(this);
1022
1023   if (field->GetMedIterator() == 0)
1024     return;
1025
1026   int ncomp = MEDfieldnComponent(FileId, field->GetMedIterator());
1027
1028   if(ncomp < 0)
1029     {
1030     field->SetNumberOfComponent(-1);
1031     vtkErrorMacro("cannot read the number of component of field "
1032         << field->GetMedIterator())
1033     return;
1034     }
1035
1036   field->SetNumberOfComponent(ncomp);
1037
1038   char* units = new char[MED_SNAME_SIZE * ncomp + 1];
1039   char* componentNames = new char[MED_SNAME_SIZE * ncomp + 1];
1040   memset(units, '\0', MED_SNAME_SIZE * ncomp + 1);
1041   memset(componentNames, '\0', MED_SNAME_SIZE * ncomp + 1);
1042
1043   //med_type_champ dataType;
1044   med_field_type dataType;
1045   med_int nstep;
1046   med_bool localmesh;
1047
1048   char name[MED_NAME_SIZE+1] = "";
1049   char meshName[MED_NAME_SIZE+1] = "";
1050   char unit[MED_SNAME_SIZE+1] = "";
1051
1052   if( MEDfieldInfo( FileId,
1053                     field->GetMedIterator(),
1054                     name,
1055                     meshName,
1056                     &localmesh,
1057                     &dataType,
1058                     componentNames,
1059                     units,
1060                     unit,
1061                     &nstep) < 0)
1062     {
1063     vtkErrorMacro("cannot read the informations on field "
1064         << field->GetMedIterator())
1065     return;
1066     }
1067
1068   field->SetName(name);
1069   field->SetMeshName(meshName);
1070   field->SetTimeUnit(unit);
1071   field->SetDataType(dataType);
1072   field->SetLocal(localmesh);
1073
1074   for(int comp = 0; comp < ncomp; comp++)
1075     {
1076     char unit[MED_NAME_SIZE + 1] = "";
1077     memcpy(unit, units + MED_SNAME_SIZE * comp, MED_SNAME_SIZE * sizeof(char));
1078     field->GetUnit()->SetValue(comp, unit);
1079
1080     char compName[MED_SNAME_SIZE + 1] = "";
1081     memcpy(compName, componentNames + MED_SNAME_SIZE * comp, MED_SNAME_SIZE
1082         * sizeof(char));
1083     field->GetComponentName()->SetValue(comp, compName);
1084     }
1085
1086   delete[] units;
1087   delete[] componentNames;
1088
1089   med_int ninterp = MEDfieldnInterp(FileId, field->GetName());
1090   if(ninterp < 0)
1091     {
1092     vtkErrorMacro("Error in MEDfieldnInterp");
1093     return;
1094     }
1095
1096   field->AllocateNumberOfInterpolation(ninterp);
1097
1098   for(med_int interpit=0; interpit<ninterp; interpit++)
1099     {
1100     vtkMedInterpolation* interp = field->GetInterpolation(interpit);
1101     interp->SetMedIterator(interpit + 1);
1102     this->ReadInterpolationInformation(interp);
1103     }
1104
1105   vtkMedFieldStep* previousStep = NULL;
1106
1107   for(med_int csit = 0; csit < nstep; csit++)
1108     {
1109     vtkMedFieldStep* step = vtkMedFieldStep::New();
1110     step->SetMedIterator(csit + 1);
1111     step->SetParentField(field);
1112     this->ReadFieldStepInformation(step, csit == 0);
1113     field->AddFieldStep(step);
1114     step->SetPreviousStep(previousStep);
1115     previousStep = step;
1116     step->Delete();
1117     }
1118 }
1119
1120 void vtkMedDriver30::ReadFieldStepInformation(vtkMedFieldStep* step, bool readAllEntityInfo)
1121 {
1122   vtkMedComputeStep cs;
1123   vtkMedComputeStep meshcs;
1124   vtkMedField* field = step->GetParentField();
1125
1126   FileOpen open(this);
1127
1128   if( MEDfieldComputingStepMeshInfo(
1129         FileId,
1130         field->GetName(),
1131         step->GetMedIterator(),
1132         &cs.TimeIt,
1133         &cs.IterationIt,
1134         &cs.TimeOrFrequency,
1135         &meshcs.TimeIt,
1136         &meshcs.IterationIt) < 0)
1137     {
1138     vtkErrorMacro("Error in MEDfieldComputingStepMeshInfo");
1139     return;
1140     }
1141
1142   step->SetComputeStep(cs);
1143   step->SetMeshComputeStep(meshcs);
1144
1145   if(!readAllEntityInfo || step->GetEntityInfoLoaded())
1146     return;
1147
1148   step->SetEntityInfoLoaded(1);
1149   
1150   vtkMedFile* file = field->GetParentFile();
1151   vtkMedMesh* mesh = file->GetMesh(field->GetMeshName());
1152   
1153   if(mesh == NULL)
1154     return;
1155   
1156   //rnv begin: fix the  "22335: [CEA 954] Paravis 7.2.0 doesn't read ELNO fields" regression.
1157   //           this piece of code needed for the reading ELNO fields
1158   std::set<vtkMedEntity> tmp_entities;
1159   std::set<vtkMedEntity> entities;
1160   mesh->GatherMedEntities(tmp_entities);
1161   
1162   std::set<vtkMedEntity>::iterator tmp_entity_it = tmp_entities.begin();
1163   while(tmp_entity_it != tmp_entities.end())
1164     {
1165       vtkMedEntity entity = *tmp_entity_it;
1166       tmp_entity_it++;
1167       entities.insert(entity);
1168       if(entity.EntityType == MED_CELL)
1169         {
1170           vtkMedEntity newEntity;
1171           newEntity.EntityType = MED_NODE_ELEMENT;
1172           newEntity.GeometryType = entity.GeometryType;
1173           newEntity.GeometryName = entity.GeometryName;
1174           entities.insert(newEntity);
1175         }
1176     }  
1177   //rnv end
1178   
1179   std::set<vtkMedEntity>::iterator entity_it = entities.begin();
1180   while(entity_it != entities.end())
1181     {
1182     vtkMedEntity entity = *entity_it;
1183     entity_it++;
1184
1185     med_int nvalues = 0;
1186     med_int nprofile;
1187     char profilename[MED_NAME_SIZE+1] = "";
1188     char localizationname[MED_NAME_SIZE+1] = "";
1189
1190     nprofile = MEDfieldnProfile(
1191         this->FileId, 
1192         field->GetName(),
1193         step->GetComputeStep().TimeIt,
1194         step->GetComputeStep().IterationIt,
1195         entity.EntityType,
1196         entity.GeometryType,
1197         profilename,
1198         localizationname);
1199     if(nprofile < 0)
1200       {
1201       vtkErrorMacro("MEDfieldnProfile");
1202       continue;
1203       }
1204
1205     bool hasprofile = (nprofile > 0);
1206     if(!hasprofile)
1207       {
1208       nprofile = 1;
1209       }
1210
1211     med_int profilesize;
1212     med_int nintegrationpoint;
1213     
1214     for(int pid=0; pid<nprofile; pid++)
1215       {
1216       med_int medid = (hasprofile ? pid+1 : -1);
1217       nvalues = MEDfieldnValueWithProfile(
1218                   this->FileId, 
1219                   field->GetName(),
1220                   step->GetComputeStep().TimeIt,
1221                   step->GetComputeStep().IterationIt,
1222                   entity.EntityType,
1223                   entity.GeometryType,
1224                   medid,
1225                   MED_COMPACT_PFLMODE,
1226                   profilename,
1227                   &profilesize,
1228                   localizationname,
1229                   &nintegrationpoint);
1230             
1231       if(nvalues < 0)
1232         {
1233         vtkErrorMacro("MEDfieldnValueWithProfile");
1234         continue;
1235         }
1236       else if(nvalues > 0)
1237         {
1238         // I have found a profile with values, stop the loop here
1239         break;
1240         }
1241       }
1242
1243     if(nvalues > 0)
1244       {
1245       vtkMedFieldOverEntity* fieldOverEntity = vtkMedFieldOverEntity::New();
1246       step->AppendFieldOverEntity(fieldOverEntity);
1247       fieldOverEntity->Delete();
1248
1249       fieldOverEntity->SetParentStep(step);
1250       fieldOverEntity->SetEntity(entity);
1251
1252       this->ReadFieldOverEntityInformation(fieldOverEntity);
1253       }
1254     }
1255 }
1256
1257 void vtkMedDriver30::ReadFieldOverEntityInformation(vtkMedFieldOverEntity* fieldOverEntity)
1258 {
1259   FileOpen open(this);
1260
1261   vtkMedFieldStep* step = fieldOverEntity->GetParentStep();
1262   vtkMedField* field = step->GetParentField();
1263   vtkMedEntity entity = fieldOverEntity->GetEntity();
1264
1265   const char* fieldName = field->GetName();
1266   const vtkMedComputeStep& cs = step->GetComputeStep();
1267
1268   char profilename[MED_NAME_SIZE+1] = "";
1269   char localizationname[MED_NAME_SIZE+1] = "";
1270
1271   med_int nProfiles = MEDfieldnProfile(
1272                         this->FileId,
1273                         fieldName,
1274                         cs.TimeIt,
1275                         cs.IterationIt,
1276                         entity.EntityType,
1277                         entity.GeometryType,
1278                         profilename,
1279                         localizationname);
1280
1281   if(nProfiles < 0)
1282     {
1283     vtkErrorMacro("MEDfieldnProfile");
1284     }
1285   else if(nProfiles == 0)
1286     {
1287     fieldOverEntity->SetHasProfile(0);
1288     nProfiles = 1;
1289     }
1290   else
1291     {
1292     fieldOverEntity->SetHasProfile(1);
1293     }
1294   fieldOverEntity->AllocateNumberOfFieldOnProfile(nProfiles);
1295   for(int profit = 0; profit < nProfiles; profit++)
1296     {
1297     vtkMedFieldOnProfile* fop = fieldOverEntity->GetFieldOnProfile(profit);
1298     med_int medid = (fieldOverEntity->GetHasProfile()? profit+1: -1);
1299     fop->SetMedIterator(medid);
1300     fop->SetParentFieldOverEntity(fieldOverEntity);
1301     this->ReadFieldOnProfileInformation(fop);
1302     }
1303 }
1304
1305 void vtkMedDriver30::ReadFieldOnProfileInformation(vtkMedFieldOnProfile* fop)
1306 {
1307   vtkMedFieldOverEntity* fieldOverEntity = fop->GetParentFieldOverEntity();
1308   vtkMedFieldStep* step = fieldOverEntity->GetParentStep();
1309   vtkMedField* field = step->GetParentField();
1310
1311   const vtkMedComputeStep& cs = step->GetComputeStep();
1312   med_int profilesize;
1313   med_int nbofintegrationpoint;
1314
1315   char profileName[MED_NAME_SIZE+1] = "";
1316   char localizationName[MED_NAME_SIZE+1] = "";
1317
1318   med_int nvalue = MEDfieldnValueWithProfile(FileId,
1319                     field->GetName(),
1320                     cs.TimeIt,
1321                     cs.IterationIt,
1322                     fieldOverEntity->GetEntity().EntityType,
1323                     fieldOverEntity->GetEntity().GeometryType,
1324                     fop->GetMedIterator(),
1325                     MED_COMPACT_STMODE,
1326                     profileName,
1327                     &profilesize,
1328                     localizationName,
1329                     &nbofintegrationpoint);
1330
1331   if(nvalue < 0)
1332     {
1333     vtkErrorMacro("Error while reading MEDfieldnValueWithProfile");
1334     }
1335
1336   fop->SetProfileName(profileName);
1337   fop->SetLocalizationName(localizationName);
1338   fop->SetNumberOfValues(nvalue);
1339   fop->SetNumberOfIntegrationPoint(nbofintegrationpoint);
1340   fop->SetProfileSize(profilesize);
1341 }
1342
1343 void vtkMedDriver30::ReadMeshInformation(vtkMedMesh* mesh)
1344 {
1345   FileOpen open(this);
1346
1347   med_int mdim = 0;
1348   med_int sdim = 0;
1349   med_mesh_type meshtype;
1350
1351   med_sorting_type sortingtype;
1352   med_int nstep = 0;
1353   med_axis_type axistype;
1354   med_int naxis;
1355
1356   if ( (naxis=MEDmeshnAxis(this->FileId, mesh->GetMedIterator())) <0 )
1357     {
1358     vtkDebugMacro("Error reading mesh axis number");
1359     }
1360
1361   if(naxis == 0)
1362     naxis=MEDmeshnAxis(this->FileId, mesh->GetMedIterator());
1363
1364   char axisname[3*MED_SNAME_SIZE+1] = "";
1365   char axisunit[3*MED_SNAME_SIZE+1] = "";
1366   char name[MED_NAME_SIZE+1] = "";
1367   char description[MED_COMMENT_SIZE+1] = "";
1368   char timeUnit[MED_SNAME_SIZE+1] = "";
1369
1370   if( MEDmeshInfo( this->FileId,
1371         mesh->GetMedIterator(),
1372         name,
1373         &sdim,
1374         &mdim,
1375         &meshtype,
1376         description,
1377         timeUnit,
1378         &sortingtype,
1379         &nstep,
1380         &axistype,
1381         axisname,
1382         axisunit ) )
1383     {
1384     vtkErrorMacro("Error reading mesh");
1385     }
1386
1387   mesh->SetName(name);
1388   mesh->SetDescription(description);
1389   mesh->SetTimeUnit(timeUnit);
1390   mesh->SetSpaceDimension(sdim);
1391   mesh->SetMeshDimension(mdim);
1392   mesh->SetMeshType(meshtype);
1393   mesh->SetSortingType(sortingtype);
1394   mesh->SetAxisType(axistype);
1395   mesh->SetNumberOfAxis(naxis);
1396
1397   for(int axis = 0; axis < naxis; axis++)
1398     {
1399     char theaxisname[MED_SNAME_SIZE+1] = "";
1400     char theaxisunit[MED_SNAME_SIZE+1] = "";
1401     strncpy(theaxisname, axisname + axis*MED_SNAME_SIZE, MED_SNAME_SIZE);
1402     strncpy(theaxisunit, axisunit + axis*MED_SNAME_SIZE, MED_SNAME_SIZE);
1403     mesh->GetAxisName()->SetValue(axis, theaxisname);
1404     mesh->GetAxisUnit()->SetValue(axis, theaxisunit);
1405     }
1406
1407   char universalName[MED_LNAME_SIZE+1] = "";
1408
1409   if(MEDmeshUniversalNameRd(this->FileId, name,
1410       universalName) < 0)
1411     {
1412     vtkDebugMacro("MEDmeshUniversalNameRd < 0");
1413     }
1414   mesh->SetUniversalName(universalName);
1415
1416   // read the Information data of all families.
1417   // writing the family 0 is optional,
1418   // but I need it, so add it if it is not here.
1419
1420   med_int nfam = MEDnFamily(this->FileId, name);
1421
1422   for(int index = 0; index < nfam; index++)
1423     {
1424     vtkMedFamily* family = vtkMedFamily::New();
1425     family->SetMedIterator(index + 1);
1426     this->ReadFamilyInformation(mesh, family);
1427     family->Delete();
1428     }
1429
1430   // this creates a family 0 if none has been read
1431   vtkMedFamily* familyZeroOnCell = mesh->GetOrCreateCellFamilyById(0);
1432   vtkMedFamily* familyZeroOnPoint = mesh->GetOrCreatePointFamilyById(0);
1433
1434   // Load Information regarding the grid type
1435   if(meshtype == MED_STRUCTURED_MESH)
1436     {
1437     // Do it for structured data
1438     med_grid_type mtg;
1439     if(MEDmeshGridTypeRd(FileId, name, &mtg) < 0)
1440       {
1441       vtkErrorMacro("Error during structured grid Information loading.");
1442       return;
1443       }
1444     mesh->SetStructuredGridType(mtg);
1445     }
1446
1447   vtkMedGrid* previousGrid = NULL;
1448   for(int gid=1; gid <= nstep; gid++)
1449     {
1450     vtkMedComputeStep cs;
1451     if(MEDmeshComputationStepInfo(FileId,
1452                                   name,
1453                                   gid,
1454                                   &cs.TimeIt,
1455                                   &cs.IterationIt,
1456                                   &cs.TimeOrFrequency) < 0)
1457       {
1458       vtkErrorMacro("MEDmeshComputationStepInfo error");
1459       }
1460     // Load Information regarding the grid type
1461     vtkMedGrid* grid = NULL;
1462     if(meshtype == MED_STRUCTURED_MESH)
1463       {
1464       switch(mesh->GetStructuredGridType())
1465         {
1466         case MED_CARTESIAN_GRID:
1467           grid = vtkMedCartesianGrid::New();
1468           break;
1469         case MED_POLAR_GRID:
1470           grid = vtkMedPolarGrid::New();
1471           break;
1472         case MED_CURVILINEAR_GRID:
1473           grid = vtkMedCurvilinearGrid::New();
1474           break;
1475         default:
1476           vtkErrorMacro("Unknown structured grid type " << mesh->GetStructuredGridType());
1477           return;
1478         }
1479       }
1480     else //(mesh->GetType() == MED_STRUCTURED_MESH)
1481       {
1482       grid = vtkMedUnstructuredGrid::New();
1483       }
1484     grid->SetParentMesh(mesh);
1485     grid->SetComputeStep(cs);
1486     this->ReadGridInformation(grid);
1487     mesh->AddGridStep(grid);
1488     grid->Delete();
1489     grid->SetPreviousGrid(previousGrid);
1490     previousGrid = grid;
1491     }
1492 }
1493
1494 void vtkMedDriver30::ReadLocalizationInformation(vtkMedLocalization* loc)
1495 {
1496   FileOpen open(this);
1497
1498   med_int ngp;
1499   med_int spaceDimension;
1500   med_geometry_type type_geo;
1501   med_geometry_type sectiongeotype;
1502   med_int nsectionmeshcell;
1503
1504   char name[MED_NAME_SIZE+1] = "";
1505   char interpolationName[MED_NAME_SIZE+1] = "";
1506   char sectionName[MED_NAME_SIZE+1] = "";
1507
1508   if(MEDlocalizationInfo(
1509       this->FileId,
1510       loc->GetMedIterator(),
1511       name,
1512       &type_geo,
1513       &spaceDimension,
1514       &ngp,
1515       interpolationName,
1516       sectionName,
1517       &nsectionmeshcell,
1518       &sectiongeotype ) < 0)
1519     {
1520     vtkErrorMacro("Reading information on quadrature points definition : "
1521         << loc->GetMedIterator());
1522     }
1523
1524   loc->SetName(name);
1525   loc->SetInterpolationName(interpolationName);
1526   loc->SetSectionName(sectionName);
1527   loc->SetNumberOfQuadraturePoint(ngp);
1528   loc->SetGeometryType(type_geo);
1529   loc->SetSpaceDimension(spaceDimension);
1530   loc->SetNumberOfCellInSection(nsectionmeshcell);
1531   loc->SetSectionGeometryType(sectiongeotype);
1532
1533   med_float *localCoordinates = new med_float[loc->GetSizeOfPointLocalCoordinates()];
1534   med_float *pqLocalCoordinates = new med_float[loc->GetSizeOfQuadraturePointLocalCoordinates()];
1535   med_float *weights = new med_float[loc->GetSizeOfWeights()];
1536
1537   if(MEDlocalizationRd(FileId,
1538       loc->GetName(),
1539       MED_FULL_INTERLACE,
1540       localCoordinates,
1541       pqLocalCoordinates,
1542       weights) < 0)
1543     {
1544     vtkErrorMacro("MEDlocalizationRd : " << loc->GetName());
1545     }
1546
1547   vtkDoubleArray* lc = loc->GetPointLocalCoordinates();
1548   vtkDoubleArray *pqlc = loc->GetQuadraturePointLocalCoordinates();
1549   vtkDoubleArray *w = loc->GetWeights();
1550
1551   lc->SetNumberOfValues(loc->GetSizeOfPointLocalCoordinates());
1552   for(int i=0; i<loc->GetSizeOfPointLocalCoordinates(); i++)
1553     {
1554     lc->SetValue(i, localCoordinates[i]);
1555     }
1556
1557   pqlc->SetNumberOfValues(loc->GetSizeOfQuadraturePointLocalCoordinates());
1558   for(int i=0; i<loc->GetSizeOfQuadraturePointLocalCoordinates(); i++)
1559     {
1560     pqlc->SetValue(i, pqLocalCoordinates[i]);
1561     }
1562
1563   w->SetNumberOfValues(loc->GetSizeOfWeights());
1564   for(int i=0; i<loc->GetSizeOfWeights(); i++)
1565     {
1566     w->SetValue(i, weights[i]);
1567     }
1568 }
1569
1570 void vtkMedDriver30::ReadInterpolationInformation(vtkMedInterpolation* interp)
1571 {
1572
1573   med_geometry_type geotype;
1574   med_bool cellnode;
1575   med_int nbofbasisfunc;
1576   med_int nbofvariable;
1577   med_int maxdegree;
1578   med_int nmaxcoef;
1579
1580   char name[MED_NAME_SIZE+1] = "";
1581
1582   if(MEDinterpInfo (this->FileId,
1583                     interp->GetMedIterator(),
1584                     name,
1585                     &geotype, &cellnode, &nbofbasisfunc,
1586                     &nbofvariable, &maxdegree, &nmaxcoef) < 0)
1587     {
1588     vtkErrorMacro("MEDinterpInfo");
1589     return;
1590     }
1591
1592   interp->SetName(name);
1593   interp->SetGeometryType(geotype);
1594   interp->SetIsCellNode(cellnode);
1595   interp->SetNumberOfVariable(nbofvariable);
1596   interp->SetMaximumDegree(maxdegree);
1597   interp->SetMaximumNumberOfCoefficient(nmaxcoef);
1598   interp->AllocateNumberOfBasisFunction(nbofbasisfunc);
1599
1600   for(int basisid=0; basisid < nbofbasisfunc; basisid++)
1601     {
1602     vtkMedFraction* func = interp->GetBasisFunction(basisid);
1603     func->SetNumberOfVariable(nbofvariable);
1604
1605     med_int ncoef = MEDinterpBaseFunctionCoefSize (
1606         this->FileId,
1607         interp->GetName(),
1608         basisid+1);
1609     func->SetNumberOfCoefficients(ncoef);
1610
1611     if(ncoef <= 0 || nbofvariable <= 0)
1612       continue;
1613
1614     med_int *power = new med_int[nbofvariable * ncoef];
1615     med_float *coefficient = new med_float[ncoef];
1616
1617     if(MEDinterpBaseFunctionRd  (
1618         this->FileId,
1619         interp->GetName(),
1620         basisid+1,
1621         &ncoef,
1622         power,
1623         coefficient) < 0)
1624       {
1625       vtkErrorMacro("MEDinterpBaseFunctionRd");
1626       continue;
1627       }
1628     vtkDoubleArray* coeffs = func->GetCoefficients();
1629     for(int cid=0; cid < ncoef; cid++)
1630       {
1631       coeffs->SetValue(cid, coefficient[cid]);
1632       }
1633     vtkIntArray* powers = func->GetPowers();
1634     for(int pid=0; pid < ncoef*nbofvariable; pid++)
1635       {
1636       powers->SetValue(pid, power[pid]);
1637       }
1638
1639     delete[] power;
1640     delete[] coefficient;
1641     }
1642 }
1643
1644 void vtkMedDriver30::ReadSupportMeshInformation(
1645     vtkMedMesh* supportMesh)
1646 {
1647   FileOpen open(this);
1648
1649   char supportmeshname[MED_NAME_SIZE+1] = "";
1650   char description[MED_COMMENT_SIZE+1] = "";
1651   med_int spacedim;
1652   med_int meshdim;
1653   med_axis_type axistype;
1654   char axisname[3*MED_SNAME_SIZE+1] = "";
1655   char axisunit[3*MED_SNAME_SIZE+1] = "";
1656
1657   if(MEDsupportMeshInfo (this->FileId,
1658                          supportMesh->GetMedIterator(),
1659                          supportmeshname,
1660                          &spacedim,
1661                          &meshdim,
1662                          description,
1663                          &axistype,
1664                          axisname,
1665                          axisunit) < 0)
1666     {
1667     vtkErrorMacro("MEDsupportMeshInfo");
1668     }
1669
1670   supportMesh->SetName(supportmeshname);
1671   supportMesh->SetDescription(description);
1672   supportMesh->SetSpaceDimension(spacedim);
1673   supportMesh->SetMeshDimension(meshdim);
1674   supportMesh->SetAxisType(axistype);
1675   for(int dim = 0; dim < 3; dim++)
1676     {
1677     char axisname_dim[MED_SNAME_SIZE+1] = "";
1678     char axisunit_dim[MED_SNAME_SIZE+1] = "";
1679
1680     strncpy(axisname_dim, axisname+dim*MED_SNAME_SIZE, MED_SNAME_SIZE);
1681     strncpy(axisunit_dim, axisunit+dim*MED_SNAME_SIZE, MED_SNAME_SIZE);
1682
1683     supportMesh->GetAxisName()->SetValue(dim, axisname_dim);
1684     supportMesh->GetAxisUnit()->SetValue(dim, axisunit_dim);
1685     }
1686
1687   return;
1688 }
1689
1690 void vtkMedDriver30::LoadFamilyIds(vtkMedEntityArray* array)
1691 {
1692   if(array->IsFamilyIdsLoaded())
1693     return;
1694
1695   FileOpen open(this);
1696
1697   vtkMedGrid* grid = array->GetParentGrid();
1698
1699   vtkMedComputeStep cs = grid->GetComputeStep();
1700
1701   // first, find if the family ids are implicit or explicit
1702   med_bool changement, transformation;
1703
1704   med_int nfamid = MEDmeshnEntity(this->FileId,
1705                       grid->GetParentMesh()->GetName(),
1706                       cs.TimeIt,
1707                       cs.IterationIt,
1708                       array->GetEntity().EntityType,
1709                       array->GetEntity().GeometryType,
1710                       MED_FAMILY_NUMBER,
1711                       MED_NO_CMODE,
1712                       &changement,
1713                       &transformation);
1714
1715   if(nfamid == array->GetNumberOfEntity())
1716     {
1717
1718     vtkMedIntArray* famIds = vtkMedIntArray::New();
1719     array->SetFamilyIds(famIds);
1720     famIds->Delete();
1721
1722     famIds->SetNumberOfTuples(nfamid);
1723
1724     if ( MEDmeshEntityFamilyNumberRd(
1725             this->FileId,
1726             grid->GetParentMesh()->GetName(),
1727             cs.TimeIt,
1728             cs.IterationIt,
1729             array->GetEntity().EntityType,
1730             array->GetEntity().GeometryType,
1731             famIds->GetPointer(0) ) < 0)
1732       {
1733       vtkWarningMacro("Error loading the family ids of entity "
1734         << array->GetEntity().EntityType
1735         << " " << array->GetEntity().GeometryType
1736         << " on mesh " << grid->GetParentMesh()->GetName());
1737       array->SetFamilyIds(NULL);
1738       }
1739     }
1740   else
1741     {
1742     vtkDebugMacro("NumberOfEntity != Number of family ids");
1743     array->SetFamilyIds(NULL);
1744     }
1745
1746   array->ComputeFamilies();
1747 }
1748
1749 void vtkMedDriver30::LoadPointGlobalIds(vtkMedGrid* grid)
1750 {
1751   if(grid->IsPointGlobalIdsLoaded())
1752     return;
1753
1754   FileOpen open(this);
1755
1756   vtkMedIntArray* globalIds = vtkMedIntArray::New();
1757   grid->SetPointGlobalIds(globalIds);
1758   globalIds->Delete();
1759
1760   globalIds->SetNumberOfTuples(grid->GetNumberOfPoints());
1761
1762   if( MEDmeshEntityNumberRd (
1763         this->FileId,
1764         grid->GetParentMesh()->GetName(),
1765         grid->GetComputeStep().TimeIt,
1766         grid->GetComputeStep().IterationIt,
1767         MED_NODE,
1768         MED_NONE,
1769         globalIds->GetPointer(0) ) < 0)
1770     {
1771     grid->SetPointGlobalIds(NULL);
1772     }
1773 }
1774
1775 void vtkMedDriver30::LoadCoordinates(vtkMedGrid* grid)
1776 {
1777   if(grid->IsCoordinatesLoaded())
1778     return;
1779
1780   vtkMedRegularGrid* rgrid = vtkMedRegularGrid::SafeDownCast(grid);
1781   if(rgrid != NULL)
1782     {
1783     this->LoadRegularGridCoordinates(rgrid);
1784     return;
1785     }
1786
1787   vtkMedUnstructuredGrid* ugrid = vtkMedUnstructuredGrid::SafeDownCast(grid);
1788   vtkMedCurvilinearGrid* cgrid = vtkMedCurvilinearGrid::SafeDownCast(grid);
1789   if(ugrid == NULL && cgrid == NULL)
1790     {
1791     //TODO : deal with structured grids
1792     vtkWarningMacro("this kind of grid is not yet supported");
1793     return;
1794     }
1795
1796   if(grid->GetUsePreviousCoordinates())
1797     {
1798     vtkMedGrid* previousgrid = grid->GetPreviousGrid();
1799     if(previousgrid == NULL)
1800       {
1801       vtkErrorMacro("coordiantes have not changed, "
1802                     << "but there is no previous grid!");
1803       return;
1804       }
1805
1806     this->LoadCoordinates(previousgrid);
1807     if(ugrid != NULL)
1808       ugrid->SetCoordinates(vtkMedUnstructuredGrid::SafeDownCast(previousgrid)
1809                             ->GetCoordinates());
1810     if(cgrid != NULL)
1811       cgrid->SetCoordinates(vtkMedCurvilinearGrid::SafeDownCast(previousgrid)
1812                             ->GetCoordinates());
1813     }
1814   else
1815     {
1816
1817     FileOpen open(this);
1818
1819     vtkDataArray* coords = vtkMedUtilities::NewCoordArray();
1820     if(ugrid != NULL)
1821       ugrid->SetCoordinates(coords);
1822     if(cgrid != NULL)
1823       cgrid->SetCoordinates(coords);
1824     coords->Delete();
1825
1826     vtkMedComputeStep cs = grid->GetComputeStep();
1827
1828     coords->SetNumberOfComponents(grid->GetParentMesh()->GetSpaceDimension());
1829     coords->SetNumberOfTuples(grid->GetNumberOfPoints());
1830
1831     if ( MEDmeshNodeCoordinateRd( this->FileId,
1832                                   grid->GetParentMesh()->GetName(),
1833                                   cs.TimeIt,
1834                                   cs.IterationIt,
1835                                   MED_FULL_INTERLACE,
1836                                   (med_float*) coords->GetVoidPointer(0) ) < 0)
1837       {
1838       vtkErrorMacro("Load Coordinates for mesh "
1839           << grid->GetParentMesh()->GetName());
1840       }
1841     }
1842 }
1843
1844 void vtkMedDriver30::LoadProfile(vtkMedProfile* profile)
1845 {
1846   if(!profile || profile->IsLoaded())
1847     return;
1848
1849   FileOpen open(this);
1850
1851   vtkMedIntArray* indices = vtkMedIntArray::New();
1852   profile->SetIds(indices);
1853   indices->Delete();
1854
1855   indices->SetNumberOfTuples(profile->GetNumberOfElement());
1856
1857   char name[MED_NAME_SIZE+1] = "";
1858
1859   if( MEDprofileRd(this->FileId,
1860                    profile->GetName(),
1861                    indices->GetPointer(0) ) < 0)
1862     {
1863     vtkErrorMacro("Reading profile indices ");
1864     }
1865 }
1866
1867 void vtkMedDriver30::LoadConnectivity(vtkMedEntityArray* array)
1868 {
1869   if(array->IsConnectivityLoaded())
1870     return;
1871
1872   FileOpen open(this);
1873
1874   vtkMedGrid* grid = array->GetParentGrid();
1875
1876   grid = array->GetParentGrid();
1877
1878   const char* meshName = grid->GetParentMesh()->GetName();
1879
1880   vtkMedIntArray* conn = vtkMedIntArray::New();
1881   array->SetConnectivityArray(conn);
1882   conn->Delete();
1883
1884   vtkMedComputeStep cs = grid->GetComputeStep();
1885
1886   med_bool change;
1887   med_bool transformation;
1888
1889   if(array->GetEntity().GeometryType == MED_POLYGON)
1890     {
1891     // first check if we have points
1892     med_int connSize = MEDmeshnEntity(
1893                             this->FileId,
1894                             meshName,
1895                             cs.TimeIt,
1896                             cs.IterationIt,
1897                             array->GetEntity().EntityType,
1898                             MED_POLYGON,
1899                             MED_CONNECTIVITY,
1900                             array->GetConnectivity(),
1901                             &change,
1902                             &transformation);
1903
1904     if (connSize < 0)
1905       {
1906       vtkErrorMacro(<< "Error while reading polygons connectivity size."
1907                     << endl );
1908             return;
1909       }
1910
1911     conn->SetNumberOfTuples(connSize);
1912
1913     // How many polygon in the mesh in nodal connectivity mode
1914     // For the polygons, we get the size of array index
1915     med_int indexsize;
1916     if ((indexsize = MEDmeshnEntity(this->FileId,
1917                                     meshName,
1918                                     cs.TimeIt,
1919                                     cs.IterationIt,
1920                                     array->GetEntity().EntityType,
1921                                     MED_POLYGON,
1922                                     MED_INDEX_NODE,
1923                                     array->GetConnectivity(),
1924                                     &change,
1925                                     &transformation )) < 0)
1926       {
1927       vtkErrorMacro(<< "Error while reading polygons array index." << endl );
1928             return;
1929       }
1930
1931     vtkMedIntArray* index = vtkMedIntArray::New();
1932     array->SetFaceIndex(index);
1933     index->Delete();
1934
1935     index->SetNumberOfTuples( indexsize );
1936
1937     if ( MEDmeshPolygonRd(this->FileId,
1938                           meshName,
1939                           cs.TimeIt,
1940                           cs.IterationIt,
1941                           array->GetEntity().EntityType,
1942                           array->GetConnectivity(),
1943                           index->GetPointer(0),
1944                           conn->GetPointer(0) ) < 0)
1945       {
1946       vtkErrorMacro(<< "MEDmeshPolygonRd");
1947       return;
1948       }
1949     }
1950   else if(array->GetEntity().GeometryType == MED_POLYHEDRON)
1951     {
1952
1953     vtkIdType connSize = MEDmeshnEntity(this->FileId,
1954                                         meshName,
1955                                         grid->GetComputeStep().TimeIt,
1956                                         grid->GetComputeStep().IterationIt,
1957                                         array->GetEntity().EntityType,
1958                                         MED_POLYHEDRON,
1959                                         MED_CONNECTIVITY,
1960                                         array->GetConnectivity(),
1961                                         &change,
1962                                         &transformation);
1963     if (connSize < 0)
1964       {
1965       vtkErrorMacro(<< "Error while reading polyhedrons connectivity size."
1966                     << endl );
1967             return;
1968       }
1969
1970     conn->SetNumberOfTuples(connSize);
1971
1972     vtkMedIntArray* faceIndex = vtkMedIntArray::New();
1973     array->SetFaceIndex(faceIndex);
1974     faceIndex->Delete();
1975
1976     vtkMedIntArray* nodeIndex = vtkMedIntArray::New();
1977     array->SetNodeIndex(nodeIndex);
1978     nodeIndex->Delete();
1979
1980     vtkIdType np = array->GetNumberOfEntity() + 1;
1981     faceIndex->SetNumberOfTuples(np);
1982
1983     med_int nodeIndexSize;
1984
1985     if ((nodeIndexSize = MEDmeshnEntity(this->FileId,
1986                                         meshName,
1987                                         grid->GetComputeStep().TimeIt,
1988                                         grid->GetComputeStep().IterationIt,
1989                                         array->GetEntity().EntityType,
1990                                         MED_POLYHEDRON,
1991                                         MED_INDEX_NODE,
1992                                         array->GetConnectivity(),
1993                                         &change,
1994                                         &transformation )) < 0)
1995       {
1996       vtkErrorMacro(<< "Error while reading polygons array index." << endl );
1997             return;
1998       }
1999
2000     nodeIndex->SetNumberOfTuples(nodeIndexSize);
2001
2002     if (MEDmeshPolyhedronRd(this->FileId,
2003                             meshName,
2004                             cs.TimeIt,
2005                             cs.IterationIt,
2006                             array->GetEntity().EntityType,
2007                             array->GetConnectivity(),
2008                             faceIndex->GetPointer(0),
2009                             nodeIndex->GetPointer(0),
2010                             conn->GetPointer(0) ) < 0)
2011       {
2012       vtkErrorMacro("Error while reading connectivity of polyhedrons");
2013       return;
2014       }
2015
2016     }
2017   else
2018     {
2019     bool doReadConnectivity = true;
2020     if(array->GetConnectivity() == MED_NODAL)
2021       {
2022       if(array->GetEntity().EntityType == MED_STRUCT_ELEMENT)
2023         {
2024         std::cout << "  -- MED_STRUCT_ELEMENT --" << std::endl;
2025         if(array->GetStructElement() == NULL)
2026           {
2027           vtkErrorMacro("Entity type = MED_STRUCT_ELEMENT, but StructElement is not set!");
2028           return;
2029           }
2030         vtkIdType ntuple = array->GetNumberOfEntity()
2031                            * array->GetStructElement()->GetConnectivitySize();
2032
2033         conn->SetNumberOfTuples(ntuple);
2034         // particles are special : connectivity is not stored in the med file
2035         if(strcmp(array->GetStructElement()->GetName(), MED_PARTICLE_NAME) == 0 )
2036           {
2037           for(vtkIdType cellId = 0; cellId < ntuple; cellId++)
2038             {
2039             conn->SetValue(cellId, cellId+1);
2040             }
2041           doReadConnectivity = false;
2042           }
2043         }
2044       else
2045         {
2046         conn->SetNumberOfTuples(array->GetNumberOfEntity()
2047             * vtkMedUtilities::GetNumberOfPoint(
2048                 array->GetEntity().GeometryType));
2049         }
2050       }
2051     else
2052       {
2053       conn->SetNumberOfTuples(array->GetNumberOfEntity()
2054           * vtkMedUtilities::GetNumberOfSubEntity(
2055               array->GetEntity().GeometryType));
2056       }
2057
2058     if  (this->ParallelFileId == -1) // also (array->GetFilter() == NULL)
2059       {
2060       if ( (MEDmeshElementConnectivityRd(
2061             this->FileId,
2062             meshName,
2063             cs.TimeIt,
2064             cs.IterationIt,
2065             array->GetEntity().EntityType,
2066             array->GetEntity().GeometryType,
2067             array->GetConnectivity(),
2068             MED_FULL_INTERLACE,
2069             conn->GetPointer(0)) ) < 0)
2070         {
2071         vtkErrorMacro("Error while load connectivity of cells "
2072             << array->GetEntity().GeometryType);
2073         }
2074       }
2075     else
2076       {
2077       med_filter filter = MED_FILTER_INIT;
2078
2079       int    start;
2080       int    stride;
2081       int    count;
2082       int    blocksize;
2083       int    lastblocksize;
2084       array->GetFilter()->GetFilterSizes(start, stride, count,
2085                                    blocksize, lastblocksize );
2086
2087       med_int nbofconstituentpervalue = vtkMedUtilities::GetNumberOfNodes(
2088                                         array->GetEntity().GeometryType);
2089
2090       if ( MEDfilterBlockOfEntityCr( this->ParallelFileId,
2091               array->GetNumberOfEntity(),
2092             1, // one is for mesh elements, more than 1 is for fields
2093               nbofconstituentpervalue,
2094             MED_ALL_CONSTITUENT,
2095             MED_FULL_INTERLACE,
2096             MED_COMPACT_STMODE,
2097             MED_NO_PROFILE,
2098             (med_size)start,
2099             (med_size)stride,
2100             (med_size)count,
2101             (med_size)blocksize,
2102             (med_size)lastblocksize,
2103             &filter ) < 0 )
2104         {
2105         vtkErrorMacro("Filter creation ");
2106         }
2107
2108         if ( (MEDmeshElementConnectivityAdvancedRd(
2109               this->ParallelFileId,
2110               meshName,
2111               cs.TimeIt,
2112               cs.IterationIt,
2113               array->GetEntity().EntityType,
2114               array->GetEntity().GeometryType,
2115               array->GetConnectivity(),
2116               &filter,
2117               conn->GetPointer(0)) ) < 0)
2118           {
2119           vtkErrorMacro("Error while load connectivity of cells "
2120               << array->GetEntity().GeometryType);
2121           }
2122
2123       if ( MEDfilterClose( &filter ) < 0)
2124           {
2125         vtkErrorMacro("ERROR : filter closing ...");
2126           }
2127
2128       }
2129     }
2130 }
2131
2132 void vtkMedDriver30::LoadCellGlobalIds(vtkMedEntityArray* array)
2133 {
2134   if(array->IsGlobalIdsLoaded())
2135     return;
2136
2137   FileOpen open(this);
2138
2139   vtkMedIntArray* globalIds = vtkMedIntArray::New();
2140   array->SetGlobalIds(globalIds);
2141   globalIds->Delete();
2142
2143   globalIds->SetNumberOfTuples(array->GetNumberOfEntity());
2144
2145   vtkMedGrid* grid = array->GetParentGrid();
2146   vtkMedComputeStep cs = grid->GetComputeStep();
2147
2148   if( MEDmeshEntityNumberRd (
2149         this->FileId,
2150         grid->GetParentMesh()->GetName(),
2151         cs.TimeIt,
2152         cs.IterationIt,
2153         array->GetEntity().EntityType,
2154         array->GetEntity().GeometryType,
2155         globalIds->GetPointer(0) ) < 0)
2156     {
2157     array->SetGlobalIds(NULL);
2158     }
2159 }
2160
2161 void vtkMedDriver30::LoadField(vtkMedFieldOnProfile* fop, med_storage_mode mode)
2162 {
2163   FileOpen open(this);
2164
2165   vtkMedFieldOverEntity* fieldOverEntity = fop->GetParentFieldOverEntity();
2166   vtkMedFieldStep *step = fieldOverEntity->GetParentStep();
2167   vtkMedField* field = step->GetParentField();
2168   const vtkMedComputeStep& cs = step->GetComputeStep();
2169
2170   vtkDataArray* data = vtkMedUtilities::NewArray(field->GetDataType());
2171   fop->SetData(data);
2172   data->Delete();
2173
2174   med_int size;
2175   if(mode == MED_COMPACT_STMODE)
2176     {
2177     size = fop->GetNumberOfValues();
2178     }
2179   else
2180     {
2181     med_int profilesize;
2182     med_int nbofintegrationpoint;
2183     char profileName[MED_NAME_SIZE+1] = "";
2184     char localizationName[MED_NAME_SIZE+1] = "";
2185     size = MEDfieldnValueWithProfile(this->FileId,
2186                 field->GetName(),
2187                 cs.TimeIt,
2188                 cs.IterationIt,
2189                 fieldOverEntity->GetEntity().EntityType,
2190                 fieldOverEntity->GetEntity().GeometryType,
2191                 fop->GetMedIterator(),
2192                 MED_GLOBAL_STMODE,
2193                 profileName,
2194                 &profilesize,
2195                 localizationName,
2196                 &nbofintegrationpoint);
2197     }
2198
2199   if(fop->GetNumberOfIntegrationPoint() > 1)
2200     {
2201     size *= fop->GetNumberOfIntegrationPoint();
2202     }
2203
2204   data->SetNumberOfComponents(field->GetNumberOfComponent());
2205   data->SetNumberOfTuples(size);
2206   if  (this->ParallelFileId == -1)
2207     {
2208     if ( MEDfieldValueWithProfileRd(
2209           this->FileId,
2210           field->GetName(),
2211           cs.TimeIt,
2212           cs.IterationIt,
2213           fieldOverEntity->GetEntity().EntityType,
2214           fieldOverEntity->GetEntity().GeometryType,
2215           mode,
2216           fop->GetProfileName(),
2217           MED_FULL_INTERLACE,
2218           MED_ALL_CONSTITUENT,
2219           (unsigned char*) data->GetVoidPointer(0) ) < 0)
2220       {
2221       vtkErrorMacro("Error on MEDfieldValueWithProfileRd");
2222       }
2223     }
2224   else
2225     {
2226   if  (field->GetFieldType() == vtkMedField::CellField)
2227     {
2228     med_filter filter = MED_FILTER_INIT;
2229
2230     int    start;
2231     int    stride;
2232     int    count;
2233     int    blocksize;
2234     int    lastblocksize;
2235     fop->GetFilter()->GetFilterSizes(start, stride, count,
2236                                  blocksize, lastblocksize );
2237
2238     if ( MEDfilterBlockOfEntityCr( this->ParallelFileId,
2239         fop->GetNumberOfValues(),
2240           1, // one is for mesh elements, more than 1 is for fields
2241           field->GetNumberOfComponent(),
2242           MED_ALL_CONSTITUENT,
2243           MED_FULL_INTERLACE,
2244           MED_COMPACT_STMODE,
2245           MED_NO_PROFILE,
2246           (med_size)start,
2247           (med_size)stride,
2248           (med_size)count,
2249           (med_size)blocksize,
2250           (med_size)lastblocksize,
2251           &filter ) < 0 )
2252       {
2253       vtkErrorMacro("Filter creation ");
2254       }
2255
2256     if ( MEDfieldValueAdvancedRd(
2257             this->ParallelFileId,
2258             field->GetName(),
2259             cs.TimeIt,
2260             cs.IterationIt,
2261             fieldOverEntity->GetEntity().EntityType,
2262             fieldOverEntity->GetEntity().GeometryType,
2263             &filter,
2264             (unsigned char*) data->GetVoidPointer(0) ) < 0)
2265         {
2266         vtkErrorMacro("Error on MEDfieldValueAdvancedRd");
2267         }
2268
2269     if ( MEDfilterClose( &filter ) < 0)
2270         {
2271       vtkErrorMacro("ERROR : filter closing ...");
2272         }
2273       }
2274   else
2275     {//TODO : option utilisateur pour desactiver ou non les champs avec profile en //
2276     if ( MEDfieldValueWithProfileRd(
2277               this->FileId,
2278               field->GetName(),
2279               cs.TimeIt,
2280               cs.IterationIt,
2281               fieldOverEntity->GetEntity().EntityType,
2282               fieldOverEntity->GetEntity().GeometryType,
2283               mode,
2284               fop->GetProfileName(),
2285               MED_FULL_INTERLACE,
2286               MED_ALL_CONSTITUENT,
2287               (unsigned char*) data->GetVoidPointer(0) ) < 0)
2288           {
2289           vtkErrorMacro("Error on MEDfieldValueWithProfileRd");
2290           }
2291     }
2292     }
2293 }
2294
2295 void vtkMedDriver30::LoadVariableAttribute(vtkMedVariableAttribute* varatt,
2296                                            vtkMedEntityArray* array)
2297 {
2298   FileOpen open(this);
2299
2300   void  *value = NULL;
2301
2302   vtkAbstractArray* valuearray = array->GetVariableAttributeValue(varatt);
2303   // first test if this is already loaded
2304   if(valuearray != NULL && valuearray->GetNumberOfTuples() > 0)
2305     return;
2306
2307   if(valuearray == NULL)
2308     {
2309     valuearray = vtkMedUtilities::NewArray(varatt->GetAttributeType());
2310     array->SetVariableAttributeValues(varatt, valuearray);
2311     valuearray->Delete();
2312     }
2313
2314   valuearray->SetNumberOfComponents(varatt->GetNumberOfComponent());
2315   valuearray->SetNumberOfTuples(array->GetNumberOfEntity());
2316   valuearray->SetName(varatt->GetName());
2317
2318   vtkSmartPointer<vtkCharArray> chararray = vtkSmartPointer<vtkCharArray>::New();
2319
2320   if(varatt->GetAttributeType() != MED_ATT_NAME)
2321     {
2322     value = valuearray->GetVoidPointer(0);
2323     }
2324   else
2325     {
2326     chararray->SetNumberOfValues(varatt->GetNumberOfComponent() *
2327                                   array->GetNumberOfEntity() *
2328                                   MED_NAME_SIZE);
2329
2330     value = chararray->GetVoidPointer(0);
2331     }
2332
2333   vtkMedComputeStep cs = array->GetParentGrid()->GetComputeStep();
2334
2335   if(MEDmeshStructElementVarAttRd(
2336       this->FileId,
2337       array->GetParentGrid()->GetParentMesh()->GetName(),
2338       cs.TimeIt,
2339       cs.IterationIt,
2340       varatt->GetParentStructElement()->GetGeometryType(),
2341       varatt->GetName(),
2342       value) < 0)
2343     {
2344
2345     if(cs.IterationIt == MED_NO_IT && cs.TimeIt == MED_NO_DT && cs.TimeOrFrequency == MED_UNDEF_DT)
2346       {
2347       vtkErrorMacro("MEDmeshStructElementVarAttRd");
2348       return;
2349       }
2350     // try to see if I can reuse
2351     // the variable attributes of the NO_DT, NO_IT compute step
2352     vtkMedComputeStep nocs;
2353     nocs.IterationIt = MED_NO_IT;
2354     nocs.TimeIt = MED_NO_DT;
2355     nocs.TimeOrFrequency = MED_UNDEF_DT;
2356     vtkMedEntityArray* nocs_array =
2357         array->GetParentGrid()->GetParentMesh()->GetGridStep(nocs)->GetEntityArray(array->GetEntity());
2358     if(nocs_array == NULL)
2359       {
2360       nocs_array = array->GetParentGrid()->GetParentMesh()->GetGridStep(0)->GetEntityArray(array->GetEntity());
2361       }
2362
2363     if(nocs_array == NULL || nocs_array == array)
2364       {
2365       // try to force load the default compute step.
2366       if(MEDmeshStructElementVarAttRd(
2367           this->FileId,
2368           array->GetParentGrid()->GetParentMesh()->GetName(),
2369           nocs.TimeIt,
2370           nocs.IterationIt,
2371           varatt->GetParentStructElement()->GetGeometryType(),
2372           varatt->GetName(),
2373           value) < 0)
2374         {
2375         vtkErrorMacro("MEDmeshStructElementVarAttRd");
2376         return;
2377         }
2378       }
2379     else
2380       {
2381       this->LoadVariableAttribute(varatt, nocs_array);
2382       array->SetVariableAttributeValues(varatt, nocs_array->GetVariableAttributeValue(varatt));
2383       return;
2384       }
2385     }
2386
2387   // If I am here, it means that I read the values
2388   if(varatt->GetAttributeType() == MED_ATT_NAME)
2389     {
2390     char current_name[MED_NAME_SIZE+1] = "";
2391     vtkStringArray* sarray = vtkStringArray::SafeDownCast(valuearray);
2392     for(vtkIdType id = 0; id < varatt->GetNumberOfComponent() *
2393                        array->GetNumberOfEntity(); id++)
2394       {
2395       memset(current_name, '\0', MED_NAME_SIZE+1);
2396       strncpy(current_name, ((char*)value) + id*MED_NAME_SIZE, MED_NAME_SIZE);
2397       sarray->SetValue(id, current_name);
2398       }
2399     }
2400
2401   return;
2402 }
2403
2404 void vtkMedDriver30::PrintSelf(ostream& os, vtkIndent indent)
2405 {
2406   this->Superclass::PrintSelf(os, indent);
2407 }