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