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