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