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