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