Salome HOME
LOT7: fix compilation at med_int==int64
[modules/paravis.git] / src / Plugins / MEDReader / IO / vtkMEDReader.cxx
1 // Copyright (C) 2010-2019  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, or (at your option) any later version.
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 // Author : Anthony Geay
20
21 #include "vtkMEDReader.h"
22 #include "vtkGenerateVectors.h"
23 #include "MEDUtilities.hxx"
24
25 #include "vtkMultiBlockDataSet.h"
26 #include "vtkInformation.h"
27 #include "vtkDataSetAttributes.h"
28 #include "vtkStringArray.h"
29 #include "vtkMutableDirectedGraph.h"
30 #include "vtkInformationStringKey.h"
31 //
32 #include "vtkUnsignedCharArray.h"
33 #include "vtkInformationVector.h"
34 #include "vtkSmartPointer.h"
35 #include "vtkVariantArray.h"
36 #include "vtkExecutive.h"
37 #include "vtkStreamingDemandDrivenPipeline.h"
38 #include "vtkMultiTimeStepAlgorithm.h"
39 #include "vtkUnstructuredGrid.h"
40 #include "vtkInformationQuadratureSchemeDefinitionVectorKey.h"
41 #include "vtkInformationDoubleVectorKey.h"
42 #include "vtkQuadratureSchemeDefinition.h"
43 #include "vtkPointData.h"
44 #include "vtkCellData.h"
45 #include "vtkCellType.h"
46 #include "vtkCellArray.h"
47 #include "vtkDoubleArray.h"
48 #include "vtkObjectFactory.h"
49 #include "vtkInformationDataObjectMetaDataKey.h"
50
51 #ifdef MEDREADER_USE_MPI
52 #include "vtkMultiProcessController.h"
53 #include "vtkPUnstructuredGridGhostCellsGenerator.h"
54 #endif
55
56 #include "MEDFileFieldRepresentationTree.hxx"
57
58 #include <map>
59 #include <string>
60 #include <vector>
61 #include <sstream>
62 #include <algorithm>
63
64 /*!
65  * This class stores properties in loading state mode (pvsm) when the MED file has not been read yet.
66  * The file is not read beacause FileName has not been informed yet ! So this class stores properties of vtkMEDReader instance that
67  * owns it and wait the vtkMEDReader::SetFileName to apply properties afterwards.
68  */
69 class PropertyKeeper
70 {
71 public:
72   PropertyKeeper(vtkMEDReader *master):IsGVActivated(false),GVValue(0),IsCMActivated(false),CMValue(0),IsGhostActivated(false),GCGCP(1),_master(master) { }
73   void assignPropertiesIfNeeded();
74   bool arePropertiesOnTreeToSetAfter() const;
75   //
76   void pushFieldStatusEntry(const char* name, int status);
77   void pushGenerateVectorsValue(int value);
78   void pushChangeModeValue(int value);
79   void pushTimesFlagsStatusEntry(const char* name, int status);
80   void pushGhost(int value);
81 protected:
82   // pool of pairs to assign in SetFieldsStatus if needed. The use case is the load using pvsm.
83   std::vector< std::pair<std::string,int> > SetFieldsStatusPairs;
84   // generate vector
85   bool IsGVActivated;
86   int GVValue;
87   // change mode
88   bool IsCMActivated;
89   int CMValue;
90   // ghost cells
91   bool IsGhostActivated;
92   int GCGCP;
93   //
94   std::vector< std::pair<std::string,int> > TimesFlagsStatusPairs;
95   vtkMEDReader *_master;
96 };
97
98 void PropertyKeeper::assignPropertiesIfNeeded()
99 {
100   if(!this->SetFieldsStatusPairs.empty())
101     {
102       for(std::vector< std::pair<std::string,int> >::const_iterator it=this->SetFieldsStatusPairs.begin();it!=this->SetFieldsStatusPairs.end();it++)
103         _master->SetFieldsStatus((*it).first.c_str(),(*it).second);
104       this->SetFieldsStatusPairs.clear();
105     }
106   if(!this->TimesFlagsStatusPairs.empty())
107     {
108       for(std::vector< std::pair<std::string,int> >::const_iterator it=this->TimesFlagsStatusPairs.begin();it!=this->TimesFlagsStatusPairs.end();it++)
109         _master->SetTimesFlagsStatus((*it).first.c_str(),(*it).second);
110       this->TimesFlagsStatusPairs.clear();
111     }
112   if(this->IsGVActivated)
113     {
114       _master->GenerateVectors(this->GVValue);
115       this->IsGVActivated=false;
116     }
117   if(this->IsCMActivated)
118     {
119       _master->ChangeMode(this->CMValue);
120       this->IsCMActivated=false;
121     }
122   if(this->IsGhostActivated)
123     {
124       _master->GhostCellGeneratorCallForPara(this->GCGCP);
125       this->IsGhostActivated=false;
126     }
127 }
128
129 void PropertyKeeper::pushFieldStatusEntry(const char* name, int status)
130 {
131   bool found(false);
132   for(std::vector< std::pair<std::string,int> >::const_iterator it=this->SetFieldsStatusPairs.begin();it!=this->SetFieldsStatusPairs.end() && !found;it++)
133     found=(*it).first==name;
134   if(!found)
135     this->SetFieldsStatusPairs.push_back(std::pair<std::string,int>(name,status));
136 }
137
138 void PropertyKeeper::pushTimesFlagsStatusEntry(const char* name, int status)
139 {
140   bool found(false);
141   for(std::vector< std::pair<std::string,int> >::const_iterator it=this->TimesFlagsStatusPairs.begin();it!=this->TimesFlagsStatusPairs.end() && !found;it++)
142     found=(*it).first==name;
143   if(!found)
144     this->TimesFlagsStatusPairs.push_back(std::pair<std::string,int>(name,status));
145 }
146
147 void PropertyKeeper::pushGenerateVectorsValue(int value)
148 {
149   this->IsGVActivated=true;
150   this->GVValue=value;
151 }
152
153 void PropertyKeeper::pushChangeModeValue(int value)
154 {
155   this->IsCMActivated=true;
156   this->CMValue=value;
157 }
158
159 void PropertyKeeper::pushGhost(int value)
160 {
161   this->IsGhostActivated=true;
162   this->GCGCP=value;
163 }
164
165 bool PropertyKeeper::arePropertiesOnTreeToSetAfter() const
166 {
167   return !SetFieldsStatusPairs.empty();
168 }
169
170 class vtkMEDReader::vtkMEDReaderInternal
171 {
172
173 public:
174   vtkMEDReaderInternal(vtkMEDReader *master):TK(0),IsMEDOrSauv(true),IsStdOrMode(false),GenerateVect(false),SIL(0),LastLev0(-1),PK(master),MyMTime(0),GCGCP(true),FirstCall0(2)
175   {
176   }
177
178   bool PluginStart0()
179   {
180     return false; // TODO Useless and buggy
181     if(FirstCall0==0)
182       return false;
183     FirstCall0--;
184     return true;
185   }
186
187   ~vtkMEDReaderInternal()
188   {
189     if(this->SIL)
190       this->SIL->Delete();
191   }
192 public:
193   MEDFileFieldRepresentationTree Tree;
194   TimeKeeper TK;
195   std::string FileName;
196   //when true the file is MED file. when false it is a Sauv file
197   bool IsMEDOrSauv;
198   //when false -> std, true -> mode. By default std (false).
199   bool IsStdOrMode;
200   //when false -> do nothing. When true cut off or extend to nbOfCompo=3 vector arrays.
201   bool GenerateVect;
202   std::string DftMeshName;
203   // Store the vtkMutableDirectedGraph that represents links between family, groups and cell types
204   vtkMutableDirectedGraph* SIL;
205   // store the lev0 id in Tree corresponding to the TIME_STEPS in the pipeline.
206   int LastLev0;
207   // The property keeper is usable only in pvsm mode.
208   PropertyKeeper PK;
209   int MyMTime;
210   std::set<std::string> _wonderful_set;// this set is used by SetFieldsStatus method to detect the fact that SetFieldsStatus has been called for all items ! Great Items are not sorted ! Why ?
211   std::map<std::string,bool> _wonderful_ref;// this map stores the state before a SetFieldsStatus status.
212   bool GCGCP;
213
214 private:
215   unsigned char FirstCall0;
216 };
217
218 vtkStandardNewMacro(vtkMEDReader);
219
220 // vtkInformationKeyMacro(vtkMEDReader, META_DATA, DataObjectMetaData); // Here we need to customize vtkMEDReader::META_DATA method
221 // start of overload of vtkInformationKeyMacro
222 static vtkInformationDataObjectMetaDataKey *vtkMEDReader_META_DATA=new vtkInformationDataObjectMetaDataKey("META_DATA","vtkMEDReader");
223
224 vtkInformationDataObjectMetaDataKey *vtkMEDReader::META_DATA()  
225 {
226   static const char ZE_KEY[]="vtkMEDReader::META_DATA";
227   vtkInformationDataObjectMetaDataKey *ret(vtkMEDReader_META_DATA);
228   MEDCoupling::GlobalDict *gd(MEDCoupling::GlobalDict::GetInstance());
229   if(!gd->hasKey(ZE_KEY))
230     {// here META_DATA is put on global var to be exchanged with other filters without dependancy of MEDReader. Please do not change ZE_KEY !
231       std::ostringstream oss; oss << ret;
232       gd->setKeyValue(ZE_KEY,oss.str());
233     }
234   return ret;
235 }
236
237 static vtkInformationGaussDoubleVectorKey *vtkMEDReader_GAUSS_DATA=new vtkInformationGaussDoubleVectorKey("GAUSS_DATA","vtkMEDReader");
238
239 vtkInformationGaussDoubleVectorKey *vtkMEDReader::GAUSS_DATA()  
240 {
241   static const char ZE_KEY[]="vtkMEDReader::GAUSS_DATA";
242   vtkInformationGaussDoubleVectorKey *ret(vtkMEDReader_GAUSS_DATA);
243   MEDCoupling::GlobalDict *gd(MEDCoupling::GlobalDict::GetInstance());
244   if(!gd->hasKey(ZE_KEY))
245     {// here META_DATA is put on global var to be exchanged with other filters without dependancy of MEDReader. Please do not change ZE_KEY !
246       vtkInformationDoubleVectorKey *ret2(ret);
247       std::ostringstream oss; oss << ret2;
248       gd->setKeyValue(ZE_KEY,oss.str());
249     }
250   return ret;
251 }
252 // end of overload of vtkInformationKeyMacro
253
254 vtkMEDReader::vtkMEDReader():Internal(new vtkMEDReaderInternal(this))
255 {
256   this->SetNumberOfInputPorts(0);
257   this->SetNumberOfOutputPorts(1);
258 }
259
260 vtkMEDReader::~vtkMEDReader()
261 {
262   delete this->Internal;
263   this->Internal = 0;
264 }
265
266 void vtkMEDReader::Reload()
267 {
268   std::string fName((const char *)this->GetFileName());
269   delete this->Internal;
270   this->Internal=new vtkMEDReaderInternal(this);
271   this->SetFileName(fName.c_str());
272 }
273
274 int vtkMEDReader::GetServerModifTime()
275 {
276   if( !this->Internal )
277     return -1;
278   return this->Internal->MyMTime;
279 }
280
281 void vtkMEDReader::GenerateVectors(int val)
282 {
283   if ( !this->Internal )
284     return;
285   
286   if(this->Internal->FileName.empty())
287     {//pvsm mode
288       this->Internal->PK.pushGenerateVectorsValue(val);
289       return ;
290     }
291   //not pvsm mode (general case)
292   bool val2((bool)val);
293   if(val2!=this->Internal->GenerateVect)
294     {
295       this->Internal->GenerateVect=val2;
296       this->Modified();
297     }
298 }
299
300 void vtkMEDReader::ChangeMode(int newMode)
301 {
302   if ( !this->Internal )
303     return;
304   
305   if(this->Internal->FileName.empty())
306     {//pvsm mode
307       this->Internal->PK.pushChangeModeValue(newMode);
308       return ;
309     }
310   //not pvsm mode (general case)
311   this->Internal->IsStdOrMode=newMode!=0;
312   this->Modified();
313 }
314
315 void vtkMEDReader::GhostCellGeneratorCallForPara(int gcgcp)
316 {
317   if ( !this->Internal )
318     return;
319   
320   if(this->Internal->FileName.empty())
321     {//pvsm mode
322       this->Internal->PK.pushGhost(gcgcp);
323       return ;
324     }
325   bool newVal(gcgcp!=0);
326   if(newVal!=this->Internal->GCGCP)
327     {
328       this->Internal->GCGCP=newVal;
329       this->Modified();
330     }
331 }
332
333 const char *vtkMEDReader::GetSeparator()
334 {
335   return MEDFileFieldRepresentationLeavesArrays::ZE_SEP;
336 }
337
338 void vtkMEDReader::SetFileName(const char *fname)
339 {
340   if(!this->Internal)
341     return;
342   try
343     {
344       this->Internal->FileName=fname;
345       std::size_t pos(this->Internal->FileName.find_last_of('.'));
346       if(pos!=std::string::npos)
347         {
348           std::string ext(this->Internal->FileName.substr(pos));
349           if(ext.find("sauv")!=std::string::npos)
350             this->Internal->IsMEDOrSauv=false;
351         }
352       if(this->Internal->Tree.getNumberOfLeavesArrays()==0)
353         {
354           int iPart(-1),nbOfParts(-1);
355 #ifdef MEDREADER_USE_MPI
356           vtkMultiProcessController *vmpc(vtkMultiProcessController::GetGlobalController());
357           if(vmpc)
358             {
359               iPart=vmpc->GetLocalProcessId();
360               nbOfParts=vmpc->GetNumberOfProcesses();
361             }
362 #endif
363           this->Internal->Tree.loadMainStructureOfFile(this->Internal->FileName.c_str(),this->Internal->IsMEDOrSauv,iPart,nbOfParts);
364           if(!this->Internal->PK.arePropertiesOnTreeToSetAfter())
365             this->Internal->Tree.activateTheFirst();//This line manually initialize the status of server (this) with the remote client.
366           this->Internal->TK.setMaxNumberOfTimeSteps(this->Internal->Tree.getMaxNumberOfTimeSteps());
367         }
368       this->Modified();
369       this->Internal->PK.assignPropertiesIfNeeded();
370     }
371   catch(INTERP_KERNEL::Exception& e)
372     {
373       delete this->Internal;
374       this->Internal=0;
375       std::ostringstream oss;
376       oss << "Exception has been thrown in vtkMEDReader::SetFileName : " << e.what() << std::endl;
377       if(this->HasObserver("ErrorEvent") )
378         this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
379       else
380         vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
381       vtkObject::BreakOnError();
382     }
383 }
384
385 char *vtkMEDReader::GetFileName()
386 {
387   if (!this->Internal)
388     return 0;
389   return const_cast<char *>(this->Internal->FileName.c_str());
390 }
391
392 int vtkMEDReader::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
393 {
394 //  std::cout << "########################################## vtkMEDReader::RequestInformation ##########################################" << std::endl;
395   if(!this->Internal)
396     return 0;
397   try
398     {
399 //      request->Print(cout);
400       vtkInformation *outInfo(outputVector->GetInformationObject(0));
401       outInfo->Set(vtkDataObject::DATA_TYPE_NAME(),"vtkMultiBlockDataSet");
402       this->UpdateSIL(request, outInfo);
403
404       // Set the meta data graph as a meta data key in the information
405       // That's all that is needed to transfer it along the pipeline
406       outInfo->Set(vtkMEDReader::META_DATA(),this->Internal->SIL);
407
408       bool dummy(false);
409       this->PublishTimeStepsIfNeeded(outInfo,dummy);
410     }
411   catch(INTERP_KERNEL::Exception& e)
412     {
413       std::ostringstream oss;
414       oss << "Exception has been thrown in vtkMEDReader::RequestInformation : " << e.what() << std::endl;
415       if(this->HasObserver("ErrorEvent") )
416         this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
417       else
418         vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
419       vtkObject::BreakOnError();
420       return 0;
421     }
422   return 1;
423 }
424
425 int vtkMEDReader::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
426 {
427 //  std::cout << "########################################## vtkMEDReader::RequestData ##########################################" << std::endl;
428   if(!this->Internal)
429     return 0;
430   try
431     {
432 //      request->Print(cout);
433       vtkInformation *outInfo(outputVector->GetInformationObject(0));
434       vtkMultiBlockDataSet *output(vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
435       bool isUpdated(false);
436       double reqTS(0.);
437       if(outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP()))
438         reqTS=outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP());
439       ExportedTinyInfo ti;
440 #ifndef MEDREADER_USE_MPI
441       this->FillMultiBlockDataSetInstance(output,reqTS,&ti);
442 #else
443       if(this->Internal->GCGCP)
444         {
445           vtkSmartPointer<vtkPUnstructuredGridGhostCellsGenerator> gcg(vtkSmartPointer<vtkPUnstructuredGridGhostCellsGenerator>::New());
446           {
447             vtkDataSet *ret(RetrieveDataSetAtTime(reqTS,&ti));
448             gcg->SetInputData(ret);
449             ret->Delete();
450           }
451           gcg->SetUseGlobalPointIds(true);
452           gcg->SetBuildIfRequired(false);
453           gcg->Update();
454           output->SetBlock(0,gcg->GetOutput());
455         }
456       else
457         this->FillMultiBlockDataSetInstance(output,reqTS,&ti);
458 #endif
459       if(!ti.empty())
460         {
461           const std::vector<double>& data(ti.getData());
462           outInfo->Set(vtkMEDReader::GAUSS_DATA(),&data[0],(int)data.size());
463           request->Append(vtkExecutive::KEYS_TO_COPY(),vtkMEDReader::GAUSS_DATA());// Thank you to SciberQuest and DIPOLE_CENTER ! Don't understand why ! In RequestInformation it does not work !
464         }
465       output->GetInformation()->Set(vtkDataObject::DATA_TIME_STEP(),reqTS);
466       // Is it really needed ? TODO
467       this->UpdateSIL(request, outInfo);
468     }
469   catch(INTERP_KERNEL::Exception& e)
470     {
471       std::cerr << "Exception has been thrown in vtkMEDReader::RequestData : " << e.what() << std::endl;
472       return 0;
473     }
474   return 1;
475 }
476
477 void vtkMEDReader::SetFieldsStatus(const char* name, int status)
478 {
479   if( !this->Internal )
480     return;
481
482   //this->Internal->_wonderful_set.insert(name);
483   if(this->Internal->FileName.empty())
484     {//pvsm mode
485       this->Internal->PK.pushFieldStatusEntry(name,status);
486       return ;
487     }
488   if(this->Internal->_wonderful_set.empty())
489     this->Internal->_wonderful_ref=this->Internal->Tree.dumpState();// start of SetFieldsStatus serie -> store ref to compare at the end of the SetFieldsStatus serie.
490   this->Internal->_wonderful_set.insert(name);
491   //not pvsm mode (general case)
492   try
493     {
494       this->Internal->Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(this->Internal->Tree.getIdHavingZeName(name),status);
495       if((int)this->Internal->_wonderful_set.size()==GetNumberOfFieldsTreeArrays())
496         {
497           if(this->Internal->_wonderful_ref!=this->Internal->Tree.dumpState())
498             {
499               if(!this->Internal->PluginStart0())
500                 {
501                   this->Modified();
502                 }
503               this->Internal->MyMTime++;
504             }
505           this->Internal->_wonderful_set.clear();
506         }
507     }
508   catch(INTERP_KERNEL::Exception& e)
509     {
510       if(!this->Internal->FileName.empty())
511         {
512           std::cerr << "vtkMEDReader::SetFieldsStatus error : " << e.what() << " *** WITH STATUS=" << status << endl;
513           return ;
514         }
515     }
516 }
517
518 int vtkMEDReader::GetNumberOfFieldsTreeArrays()
519 {
520   if(!this->Internal)
521     return 0;
522   return this->Internal->Tree.getNumberOfLeavesArrays();
523 }
524
525 const char *vtkMEDReader::GetFieldsTreeArrayName(int index)
526 {
527   if(!this->Internal)
528     return 0;
529   return this->Internal->Tree.getNameOfC(index);
530 }
531
532 int vtkMEDReader::GetFieldsTreeArrayStatus(const char *name)
533 {
534   if(!this->Internal)
535     return -1;
536
537   int zeId(this->Internal->Tree.getIdHavingZeName(name));
538   int ret(this->Internal->Tree.getStatusOf(zeId));
539   return ret;
540 }
541
542 void vtkMEDReader::SetTimesFlagsStatus(const char *name, int status)
543 {
544   if (!this->Internal)
545     return;
546   
547   if(this->Internal->FileName.empty())
548     {//pvsm mode
549       this->Internal->PK.pushTimesFlagsStatusEntry(name,status);
550       return ;
551     }
552   //not pvsm mode (general case)
553   int pos(0);
554   std::istringstream iss(name); iss >> pos;
555   this->Internal->TK.getTimesFlagArray()[pos].first=(bool)status;
556   if(pos==(int)this->Internal->TK.getTimesFlagArray().size()-1)
557     if(!this->Internal->PluginStart0())
558       {
559         this->Modified();
560         //this->Internal->TK.printSelf(std::cerr);
561       }
562 }
563
564 int vtkMEDReader::GetNumberOfTimesFlagsArrays()
565 {
566   if(!this->Internal)
567     return 0;
568   return (int)this->Internal->TK.getTimesFlagArray().size();
569 }
570
571 const char *vtkMEDReader::GetTimesFlagsArrayName(int index)
572 {
573   return this->Internal->TK.getTimesFlagArray()[index].second.c_str();
574 }
575
576 int vtkMEDReader::GetTimesFlagsArrayStatus(const char *name)
577 {
578   if(!this->Internal)
579     return -1;
580   int pos(0);
581   std::istringstream iss(name); iss >> pos;
582   return (int)this->Internal->TK.getTimesFlagArray()[pos].first;
583 }
584
585 void vtkMEDReader::UpdateSIL(vtkInformation* request, vtkInformation *info)
586 {
587   if(!this->Internal)
588       return;
589   std::string meshName(this->Internal->Tree.getActiveMeshName());
590   if(!this->Internal->SIL || meshName!=this->Internal->DftMeshName)
591     {
592       vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::New());
593       this->BuildSIL(sil);
594       if(this->Internal->SIL)
595         this->Internal->SIL->Delete();
596       this->Internal->SIL=sil;
597       this->Internal->DftMeshName=meshName;
598     }
599 }
600
601 /*!
602  * The returned string is the name of the mesh activated which groups and families are in \a sil.
603  */
604 std::string vtkMEDReader::BuildSIL(vtkMutableDirectedGraph* sil)
605 {
606   if (!this->Internal)
607     return std::string();
608   sil->Initialize();
609   vtkSmartPointer<vtkVariantArray> childEdge(vtkSmartPointer<vtkVariantArray>::New());
610   childEdge->InsertNextValue(0);
611   vtkSmartPointer<vtkVariantArray> crossEdge(vtkSmartPointer<vtkVariantArray>::New());
612   crossEdge->InsertNextValue(1);
613   // CrossEdge is an edge linking hierarchies.
614   vtkUnsignedCharArray* crossEdgesArray=vtkUnsignedCharArray::New();
615   crossEdgesArray->SetName("CrossEdges");
616   sil->GetEdgeData()->AddArray(crossEdgesArray);
617   crossEdgesArray->Delete();
618   std::vector<std::string> names;
619   // Now build the hierarchy.
620   vtkIdType rootId=sil->AddVertex();
621   names.push_back("SIL");
622   // Add global fields root
623   vtkIdType fieldsRoot(sil->AddChild(rootId,childEdge));
624   names.push_back("FieldsStatusTree");
625   this->Internal->Tree.feedSIL(sil,fieldsRoot,childEdge,names);
626   vtkIdType meshesFamsGrpsRoot(sil->AddChild(rootId,childEdge));
627   names.push_back("MeshesFamsGrps");
628   std::string dftMeshName(this->Internal->Tree.feedSILForFamsAndGrps(sil,meshesFamsGrpsRoot,childEdge,names));
629   // This array is used to assign names to nodes.
630   vtkStringArray *namesArray(vtkStringArray::New());
631   namesArray->SetName("Names");
632   namesArray->SetNumberOfTuples(sil->GetNumberOfVertices());
633   sil->GetVertexData()->AddArray(namesArray);
634   namesArray->Delete();
635   std::vector<std::string>::const_iterator iter;
636   vtkIdType cc;
637   for(cc=0, iter=names.begin(); iter!=names.end(); ++iter, ++cc)
638     namesArray->SetValue(cc,(*iter).c_str());
639   return dftMeshName;
640 }
641
642 double vtkMEDReader::PublishTimeStepsIfNeeded(vtkInformation *outInfo, bool& isUpdated)
643 {
644   if(!this->Internal)
645     return 0.0;
646
647   int lev0(-1);
648   std::vector<double> tsteps;
649   if(!this->Internal->IsStdOrMode)
650     tsteps=this->Internal->Tree.getTimeSteps(lev0,this->Internal->TK);
651   else
652     { tsteps.resize(1); tsteps[0]=0.; }
653   isUpdated=false;
654   if(lev0!=this->Internal->LastLev0)
655     {
656       isUpdated=true;
657       double timeRange[2];
658       timeRange[0]=tsteps.front();
659       timeRange[1]=tsteps.back();
660       outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),&tsteps[0],(int)tsteps.size());
661       outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(),timeRange,2);
662       this->Internal->LastLev0=lev0;
663     }
664   return tsteps.front();
665 }
666
667 void vtkMEDReader::FillMultiBlockDataSetInstance(vtkMultiBlockDataSet *output, double reqTS, ExportedTinyInfo *internalInfo)
668 {
669   if( !this->Internal )
670     return;
671   vtkDataSet *ret(RetrieveDataSetAtTime(reqTS,internalInfo));
672   output->SetBlock(0,ret);
673   ret->Delete();
674 }
675
676 vtkDataSet *vtkMEDReader::RetrieveDataSetAtTime(double reqTS, ExportedTinyInfo *internalInfo)
677 {
678   if( !this->Internal )
679     return 0;
680   std::string meshName;
681   vtkDataSet *ret(this->Internal->Tree.buildVTKInstance(this->Internal->IsStdOrMode,reqTS,meshName,this->Internal->TK,internalInfo));
682   if(this->Internal->GenerateVect)
683     {
684       vtkGenerateVectors::Operate(ret->GetPointData());
685       vtkGenerateVectors::Operate(ret->GetCellData());
686       vtkGenerateVectors::Operate(ret->GetFieldData());
687       // The operations above have potentially created new arrays -> This breaks the optimization of StaticMesh that expects the same field arrays over time.
688       // To enforce the cache recomputation declare modification of mesh.
689       //vtkGenerateVectors::ChangeMeshTimeToUpdateCache(ret);
690     }
691   return ret;
692 }
693
694 void vtkMEDReader::PrintSelf(ostream& os, vtkIndent indent)
695 {
696   this->Superclass::PrintSelf(os, indent);
697 }