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