]> SALOME platform Git repositories - modules/paravis.git/blob - src/Plugins/MEDReader/IO/vtkMEDReader.cxx
Salome HOME
Name change of the Integration points filter menu.
[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 #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)
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 private:
192   unsigned char FirstCall0;
193 };
194
195 vtkStandardNewMacro(vtkMEDReader);
196
197 vtkMEDReader::vtkMEDReader():Internal(new vtkMEDReaderInternal(this))
198 {
199   this->SetNumberOfInputPorts(0);
200   this->SetNumberOfOutputPorts(1);
201 }
202
203 vtkMEDReader::~vtkMEDReader()
204 {
205   delete this->Internal;
206 }
207
208 void vtkMEDReader::Reload(int a)
209 {
210   if(a==0)
211     return;
212   std::cerr << "vtkMEDReader::Reload" << a << std::endl;
213   std::string fName((const char *)this->GetFileName());
214   delete this->Internal;
215   this->Internal=new vtkMEDReaderInternal(this);
216   this->SetFileName(fName.c_str());
217 }
218
219 int vtkMEDReader::GetServerModifTime()
220 {
221   return this->Internal->MyMTime;
222 }
223
224 void vtkMEDReader::GenerateVectors(int val)
225 {
226   if(this->Internal->FileName.empty())
227     {//pvsm mode
228       this->Internal->PK.pushGenerateVectorsValue(val);
229       return ;
230     }
231   //not pvsm mode (general case)
232   bool val2((bool)val);
233   if(val2!=this->Internal->GenerateVect)
234     {
235       this->Internal->GenerateVect=val2;
236       this->Modified();
237     }
238 }
239
240 void vtkMEDReader::ChangeMode(int newMode)
241 {
242   if(this->Internal->FileName.empty())
243     {//pvsm mode
244       this->Internal->PK.pushChangeModeValue(newMode);
245       return ;
246     }
247   //not pvsm mode (general case)
248   this->Internal->IsStdOrMode=newMode!=0;
249   //std::cerr << "vtkMEDReader::ChangeMode : " << this->Internal->IsStdOrMode << std::endl;
250   this->Modified();
251 }
252
253 const char *vtkMEDReader::GetSeparator()
254 {
255   return MEDFileFieldRepresentationLeavesArrays::ZE_SEP;
256 }
257
258 void vtkMEDReader::SetFileName(const char *fname)
259 {
260   try
261     {
262       this->Internal->FileName=fname;
263       std::size_t pos(this->Internal->FileName.find_last_of('.'));
264       if(pos!=std::string::npos)
265         {
266           std::string ext(this->Internal->FileName.substr(pos));
267           if(ext.find("sauv")!=std::string::npos)
268             this->Internal->IsMEDOrSauv=false;
269         }
270       if(this->Internal->Tree.getNumberOfLeavesArrays()==0)
271         {
272           int iPart(-1),nbOfParts(-1);
273 #ifdef MEDREADER_USE_MPI
274           vtkMultiProcessController *vmpc(vtkMultiProcessController::GetGlobalController());
275           if(vmpc)
276             {
277               iPart=vmpc->GetLocalProcessId();
278               nbOfParts=vmpc->GetNumberOfProcesses();
279             }
280 #endif
281           this->Internal->Tree.loadMainStructureOfFile(this->Internal->FileName.c_str(),this->Internal->IsMEDOrSauv,iPart,nbOfParts);
282           if(!this->Internal->PK.arePropertiesOnTreeToSetAfter())
283             this->Internal->Tree.activateTheFirst();//This line manually initialize the status of server (this) with the remote client.
284           this->Internal->TK.setMaxNumberOfTimeSteps(this->Internal->Tree.getMaxNumberOfTimeSteps());
285         }
286       this->Modified();
287       this->Internal->PK.assignPropertiesIfNeeded();
288     }
289   catch(INTERP_KERNEL::Exception& e)
290     {
291       delete this->Internal;
292       this->Internal=0;
293       std::ostringstream oss;
294       oss << "Exception has been thrown in vtkMEDReader::SetFileName : " << e.what() << std::endl;
295       if(this->HasObserver("ErrorEvent") )
296         this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
297       else
298         vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
299       vtkObject::BreakOnError();
300     }
301 }
302
303 char *vtkMEDReader::GetFileName()
304 {
305   return const_cast<char *>(this->Internal->FileName.c_str());
306 }
307
308 int vtkMEDReader::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
309 {
310   //std::cerr << "########################################## RequestInformation ##########################################" << std::endl;
311   if(!this->Internal)
312     return 0;
313   try
314     {
315       vtkInformation *outInfo(outputVector->GetInformationObject(0));
316       outInfo->Set(vtkDataObject::DATA_TYPE_NAME(),"vtkMultiBlockDataSet");
317       this->UpdateSIL(outInfo);
318       //
319       bool dummy(false);
320       this->PublishTimeStepsIfNeeded(outInfo,dummy);
321     }
322   catch(INTERP_KERNEL::Exception& e)
323     {
324       std::ostringstream oss;
325       oss << "Exception has been thrown in vtkMEDReader::RequestInformation : " << e.what() << std::endl;
326       if(this->HasObserver("ErrorEvent") )
327         this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
328       else
329         vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
330       vtkObject::BreakOnError();
331       return 0;
332     }
333   return 1;
334 }
335
336 int vtkMEDReader::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
337 {
338   //std::cerr << "########################################## RequestData        ##########################################";
339   if(!this->Internal)
340     return 0;
341   try
342     {
343       vtkInformation *outInfo(outputVector->GetInformationObject(0));
344       vtkMultiBlockDataSet *output(vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
345       bool isUpdated(false);
346       double reqTS(0.);
347       if(outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP()))
348         reqTS=outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP());
349       //std::cerr << reqTS << std::endl;
350       this->FillMultiBlockDataSetInstance(output,reqTS);
351       output->GetInformation()->Set(vtkDataObject::DATA_TIME_STEP(),reqTS);
352       this->UpdateSIL(outInfo);
353       //this->UpdateProgress((float) progress/((float) maxprogress-1));
354     }
355   catch(INTERP_KERNEL::Exception& e)
356     {
357       std::cerr << "Exception has been thrown in vtkMEDReader::RequestData : " << e.what() << std::endl;
358       return 0;
359     }
360   return 1;
361 }
362
363 void vtkMEDReader::SetFieldsStatus(const char* name, int status)
364 {
365   //this->Internal->_wonderful_set.insert(name);
366   if(this->Internal->FileName.empty())
367     {//pvsm mode
368       this->Internal->PK.pushFieldStatusEntry(name,status);
369       return ;
370     }
371   if(this->Internal->_wonderful_set.empty())
372     this->Internal->_wonderful_ref=this->Internal->Tree.dumpState();// start of SetFieldsStatus serie -> store ref to compare at the end of the SetFieldsStatus serie.
373   this->Internal->_wonderful_set.insert(name);
374   //not pvsm mode (general case)
375   try
376     {
377       this->Internal->Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(this->Internal->Tree.getIdHavingZeName(name),status);
378       if(this->Internal->_wonderful_set.size()==GetNumberOfFieldsTreeArrays())
379         {
380           if(this->Internal->_wonderful_ref!=this->Internal->Tree.dumpState())
381             {
382               if(!this->Internal->PluginStart0())
383                 {
384                   this->Modified();
385                 }
386               this->Internal->MyMTime++;
387             }
388           this->Internal->_wonderful_set.clear();
389         }
390     }
391   catch(INTERP_KERNEL::Exception& e)
392     {
393       if(!this->Internal->FileName.empty())
394         {
395           std::cerr << "vtkMEDReader::SetFieldsStatus error : " << e.what() << " *** WITH STATUS=" << status << endl;
396           return ;
397         }
398     }
399 }
400
401 int vtkMEDReader::GetNumberOfFieldsTreeArrays()
402 {
403   if(!this->Internal)
404     return 0;
405   return this->Internal->Tree.getNumberOfLeavesArrays();
406   //std::cerr << "vtkMEDReader::GetNumberOfFieldsTreeArrays called ! " << ret << std::endl;
407 }
408
409 const char *vtkMEDReader::GetFieldsTreeArrayName(int index)
410 {
411   return this->Internal->Tree.getNameOfC(index);
412   //std::cerr << "vtkMEDReader::GetFieldsTreeArrayName(" << index << ") called ! " << ret << std::endl;
413 }
414
415 int vtkMEDReader::GetFieldsTreeArrayStatus(const char *name)
416 {
417   int zeId(this->Internal->Tree.getIdHavingZeName(name));
418   int ret(this->Internal->Tree.getStatusOf(zeId));
419   return ret;
420 }
421
422 void vtkMEDReader::SetTimesFlagsStatus(const char *name, int status)
423 {
424   if(this->Internal->FileName.empty())
425     {//pvsm mode
426       this->Internal->PK.pushTimesFlagsStatusEntry(name,status);
427       return ;
428     }
429   //not pvsm mode (general case)
430   int pos(0);
431   std::istringstream iss(name); iss >> pos;
432   this->Internal->TK.getTimesFlagArray()[pos].first=(bool)status;
433   if(pos==this->Internal->TK.getTimesFlagArray().size()-1)
434     if(!this->Internal->PluginStart0())
435       {
436         this->Modified();
437         //this->Internal->TK.printSelf(std::cerr);
438       }
439 }
440
441 int vtkMEDReader::GetNumberOfTimesFlagsArrays()
442 {
443   if(!this->Internal)
444     return 0;
445   return (int)this->Internal->TK.getTimesFlagArray().size();
446 }
447
448 const char *vtkMEDReader::GetTimesFlagsArrayName(int index)
449 {
450   return this->Internal->TK.getTimesFlagArray()[index].second.c_str();
451 }
452
453 int vtkMEDReader::GetTimesFlagsArrayStatus(const char *name)
454 {
455   int pos(0);
456   std::istringstream iss(name); iss >> pos;
457   return (int)this->Internal->TK.getTimesFlagArray()[pos].first;
458 }
459
460 void vtkMEDReader::UpdateSIL(vtkInformation *info)
461 {
462   if(!this->Internal)
463       return ;
464   vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::New());
465   std::string meshName(this->BuildSIL(sil));
466   if(meshName!=this->Internal->DftMeshName)
467     {
468       if(this->Internal->SIL)
469         this->Internal->SIL->Delete();
470       this->Internal->SIL=sil;
471       this->Internal->DftMeshName=meshName;
472       info->Set(vtkDataObject::SIL(),this->Internal->SIL);
473       //request->AppendUnique(vtkExecutive::KEYS_TO_COPY(),vtkDataObject::SIL());
474     }
475   else
476     {
477       sil->Delete();
478     }
479 }
480
481 /*!
482  * The returned string is the name of the mesh activated which groups and families are in \a sil.
483  */
484 std::string vtkMEDReader::BuildSIL(vtkMutableDirectedGraph* sil)
485 {
486   sil->Initialize();
487   vtkSmartPointer<vtkVariantArray> childEdge(vtkSmartPointer<vtkVariantArray>::New());
488   childEdge->InsertNextValue(0);
489   vtkSmartPointer<vtkVariantArray> crossEdge(vtkSmartPointer<vtkVariantArray>::New());
490   crossEdge->InsertNextValue(1);
491   // CrossEdge is an edge linking hierarchies.
492   vtkUnsignedCharArray* crossEdgesArray=vtkUnsignedCharArray::New();
493   crossEdgesArray->SetName("CrossEdges");
494   sil->GetEdgeData()->AddArray(crossEdgesArray);
495   crossEdgesArray->Delete();
496   std::vector<std::string> names;
497   // Now build the hierarchy.
498   vtkIdType rootId=sil->AddVertex();
499   names.push_back("SIL");
500   // Add global fields root
501   vtkIdType fieldsRoot(sil->AddChild(rootId,childEdge));
502   names.push_back("FieldsStatusTree");
503   this->Internal->Tree.feedSIL(sil,fieldsRoot,childEdge,names);
504   vtkIdType meshesFamsGrpsRoot(sil->AddChild(rootId,childEdge));
505   names.push_back("MeshesFamsGrps");
506   std::string dftMeshName(this->Internal->Tree.feedSILForFamsAndGrps(sil,meshesFamsGrpsRoot,childEdge,names));
507   // This array is used to assign names to nodes.
508   vtkStringArray *namesArray(vtkStringArray::New());
509   namesArray->SetName("Names");
510   namesArray->SetNumberOfTuples(sil->GetNumberOfVertices());
511   sil->GetVertexData()->AddArray(namesArray);
512   namesArray->Delete();
513   std::vector<std::string>::const_iterator iter;
514   vtkIdType cc;
515   for(cc=0, iter=names.begin(); iter!=names.end(); ++iter, ++cc)
516     namesArray->SetValue(cc,(*iter).c_str());
517   return dftMeshName;
518 }
519
520 double vtkMEDReader::PublishTimeStepsIfNeeded(vtkInformation *outInfo, bool& isUpdated)
521 {
522   int lev0(-1);
523   std::vector<double> tsteps;
524   if(!this->Internal->IsStdOrMode)
525     tsteps=this->Internal->Tree.getTimeSteps(lev0,this->Internal->TK);
526   else
527     { tsteps.resize(1); tsteps[0]=0.; }
528   isUpdated=false;
529   if(lev0!=this->Internal->LastLev0)
530     {
531       isUpdated=true;
532       double timeRange[2];
533       timeRange[0]=tsteps.front();
534       timeRange[1]=tsteps.back();
535       outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),&tsteps[0],tsteps.size());
536       outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(),timeRange,2);
537       this->Internal->LastLev0=lev0;
538     }
539   return tsteps.front();
540 }
541
542 void vtkMEDReader::FillMultiBlockDataSetInstance(vtkMultiBlockDataSet *output, double reqTS)
543 {
544   std::string meshName;
545   vtkDataSet *ret(this->Internal->Tree.buildVTKInstance(this->Internal->IsStdOrMode,reqTS,meshName,this->Internal->TK));
546   if(this->Internal->GenerateVect)
547     {
548       vtkGenerateVectors::Operate(ret->GetPointData());
549       vtkGenerateVectors::Operate(ret->GetCellData());
550       vtkGenerateVectors::Operate(ret->GetFieldData());
551     }
552   output->SetBlock(0,ret);
553   ret->Delete();
554 }
555
556 void vtkMEDReader::PrintSelf(ostream& os, vtkIndent indent)
557 {
558   this->Superclass::PrintSelf(os, indent);
559 }