Salome HOME
new MEDReader
[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.
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 class vtkMEDReader::vtkMEDReaderInternal
56 {
57
58 public:
59   vtkMEDReaderInternal():TK(0),IsMEDOrSauv(true),IsStdOrMode(false),GenerateVect(false),SIL(0),LastLev0(-1),FirstCall0(2)
60   {
61   }
62   
63   bool PluginStart0()
64   {
65     if(FirstCall0==0)
66       return false;
67     FirstCall0--;
68     return true;
69   }
70   
71   ~vtkMEDReaderInternal()
72   {
73     if(this->SIL)
74       this->SIL->Delete();
75   }
76 public:
77   MEDFileFieldRepresentationTree Tree;
78   TimeKeeper TK;
79   std::string FileName;
80   //when true the file is MED file. when false it is a Sauv file
81   bool IsMEDOrSauv;
82   //when false -> std, true -> mode. By default std (false).
83   bool IsStdOrMode;
84   //when false -> do nothing. When true cut off or extend to nbOfCompo=3 vector arrays.
85   bool GenerateVect;
86   std::string DftMeshName;
87   // Store the vtkMutableDirectedGraph that represents links between family, groups and cell types
88   vtkMutableDirectedGraph* SIL;
89   // store the lev0 id in Tree corresponding to the TIME_STEPS in the pipeline.
90   int LastLev0;
91 private:
92   unsigned char FirstCall0;
93 };
94
95 vtkStandardNewMacro(vtkMEDReader);
96
97 vtkMEDReader::vtkMEDReader():Internal(new vtkMEDReaderInternal)
98 {
99   this->SetNumberOfInputPorts(0);
100   this->SetNumberOfOutputPorts(1);
101 }
102
103 vtkMEDReader::~vtkMEDReader()
104 {
105   delete this->Internal;
106 }
107
108 void vtkMEDReader::Reload(int a)
109 {
110   if(a==0)
111     return;
112   std::cerr << "vtkMEDReader::Reload" << a << std::endl;
113   std::string fName((const char *)this->GetFileName());
114   delete this->Internal;
115   this->Internal=new vtkMEDReaderInternal;
116   this->SetFileName(fName.c_str());
117 }
118
119 void vtkMEDReader::GenerateVectors(int val)
120 {
121   bool val2((bool)val);
122   if(val2!=this->Internal->GenerateVect)
123     {
124       this->Internal->GenerateVect=val2;
125       this->Modified();
126     }
127 }
128
129 void vtkMEDReader::ChangeMode(int newMode)
130 {
131   this->Internal->IsStdOrMode=newMode!=0;
132   //std::cerr << "vtkMEDReader::ChangeMode : " << this->Internal->IsStdOrMode << std::endl;
133   this->Modified();
134 }
135
136 const char *vtkMEDReader::GetSeparator()
137 {
138   return MEDFileFieldRepresentationLeavesArrays::ZE_SEP;
139 }
140
141 void vtkMEDReader::SetFileName(const char *fname)
142 {
143   this->Internal->FileName=fname;
144   std::size_t pos(this->Internal->FileName.find_last_of('.'));
145   if(pos!=std::string::npos)
146     {
147       std::string ext(this->Internal->FileName.substr(pos));
148       if(ext.find("sauv")!=std::string::npos)
149         this->Internal->IsMEDOrSauv=false;
150     }
151   if(this->Internal->Tree.getNumberOfLeavesArrays()==0)
152     {
153       this->Internal->Tree.loadMainStructureOfFile(this->Internal->FileName.c_str(),this->Internal->IsMEDOrSauv);
154       this->Internal->Tree.activateTheFirst();//This line manually initialize the status of server (this) with the remote client.
155       this->Internal->TK.setMaxNumberOfTimeSteps(this->Internal->Tree.getMaxNumberOfTimeSteps());
156     }
157   this->Modified();
158 }
159
160 char *vtkMEDReader::GetFileName()
161 {
162   return const_cast<char *>(this->Internal->FileName.c_str());
163 }
164
165 int vtkMEDReader::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
166 {
167   //std::cerr << "########################################## RequestInformation ##########################################" << std::endl;
168   try
169     {
170       vtkInformation *outInfo(outputVector->GetInformationObject(0));
171       outInfo->Set(vtkStreamingDemandDrivenPipeline::MAXIMUM_NUMBER_OF_PIECES(),-1);
172       outInfo->Set(vtkDataObject::DATA_TYPE_NAME(),"vtkMultiBlockDataSet");
173       this->UpdateSIL(outInfo);
174       //
175       bool dummy(false);
176       this->PublishTimeStepsIfNeeded(outInfo,dummy);
177     }
178   catch(INTERP_KERNEL::Exception& e)
179     {
180       std::cerr << "Exception has been thrown in vtkMEDReader::RequestInformation : " << e.what() << std::endl;
181       return 0;
182     }
183   return 1;
184 }
185
186 int vtkMEDReader::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
187 {
188   //std::cerr << "########################################## RequestData        ##########################################";
189   try
190     {
191       vtkInformation *outInfo(outputVector->GetInformationObject(0));
192       vtkMultiBlockDataSet *output(vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
193       bool isUpdated(false);
194       double reqTS(0.);
195       if(outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP()))
196         reqTS=outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP());
197       //std::cerr << reqTS << std::endl;
198       this->FillMultiBlockDataSetInstance(output,reqTS);
199       output->GetInformation()->Set(vtkDataObject::DATA_TIME_STEP(),reqTS);
200       this->UpdateSIL(outInfo);
201       //this->UpdateProgress((float) progress/((float) maxprogress-1));
202     }
203   catch(INTERP_KERNEL::Exception& e)
204     {
205       std::cerr << "Exception has been thrown in vtkMEDReader::RequestInformation : " << e.what() << std::endl;
206       return 0;
207     }
208   return 1;
209 }
210
211 void vtkMEDReader::SetFieldsStatus(const char* name, int status)
212 {
213   //std::cerr << "vtkMEDReader::SetFieldsStatus(" << name << "," << status << ") called !" << std::endl;
214   this->Internal->Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(this->Internal->Tree.getIdHavingZeName(name),status);
215   if(std::string(name)==GetFieldsTreeArrayName(GetNumberOfFieldsTreeArrays()-1))
216     if(!this->Internal->PluginStart0())
217       this->Modified();
218 }
219
220 int vtkMEDReader::GetNumberOfFieldsTreeArrays()
221 {
222   int ret(this->Internal->Tree.getNumberOfLeavesArrays());
223   //std::cerr << "vtkMEDReader::GetNumberOfFieldsTreeArrays called ! " << ret << std::endl;
224   return ret;
225 }
226
227 const char *vtkMEDReader::GetFieldsTreeArrayName(int index)
228 {
229   std::string ret(this->Internal->Tree.getNameOf(index));
230   //std::cerr << "vtkMEDReader::GetFieldsTreeArrayName(" << index << ") called ! " << ret << std::endl;
231   return ret.c_str();
232 }
233
234 int vtkMEDReader::GetFieldsTreeArrayStatus(const char *name)
235 {
236   int zeId(this->Internal->Tree.getIdHavingZeName(name));
237   int ret(this->Internal->Tree.getStatusOf(zeId));
238   return ret;
239 }
240
241 void vtkMEDReader::SetTimesFlagsStatus(const char *name, int status)
242 {
243   int pos(0);
244   std::istringstream iss(name); iss >> pos;
245   this->Internal->TK.getTimesFlagArray()[pos].first=(bool)status;
246   if(pos==this->Internal->TK.getTimesFlagArray().size()-1)
247     if(!this->Internal->PluginStart0())
248       {
249         this->Modified();
250         //this->Internal->TK.printSelf(std::cerr);
251       }
252 }
253
254 int vtkMEDReader::GetNumberOfTimesFlagsArrays()
255 {
256   return (int)this->Internal->TK.getTimesFlagArray().size();
257 }
258
259 const char *vtkMEDReader::GetTimesFlagsArrayName(int index)
260 {
261   return this->Internal->TK.getTimesFlagArray()[index].second.c_str();
262 }
263
264 int vtkMEDReader::GetTimesFlagsArrayStatus(const char *name)
265 {
266   int pos(0);
267   std::istringstream iss(name); iss >> pos;
268   return (int)this->Internal->TK.getTimesFlagArray()[pos].first;
269 }
270
271 void vtkMEDReader::UpdateSIL(vtkInformation *info)
272 {
273   vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::New());
274   std::string meshName(this->BuildSIL(sil));
275   if(meshName!=this->Internal->DftMeshName)
276     {
277       if(this->Internal->SIL)
278         this->Internal->SIL->Delete();
279       this->Internal->SIL=sil;
280       this->Internal->DftMeshName=meshName;
281       info->Set(vtkDataObject::SIL(),this->Internal->SIL);
282       //request->AppendUnique(vtkExecutive::KEYS_TO_COPY(),vtkDataObject::SIL());
283     }
284   else
285     {
286       sil->Delete();
287     }
288 }
289
290 /*!
291  * The returned string is the name of the mesh activated which groups and families are in \a sil.
292  */
293 std::string vtkMEDReader::BuildSIL(vtkMutableDirectedGraph* sil)
294 {
295   sil->Initialize();
296   vtkSmartPointer<vtkVariantArray> childEdge(vtkSmartPointer<vtkVariantArray>::New());
297   childEdge->InsertNextValue(0);
298   vtkSmartPointer<vtkVariantArray> crossEdge(vtkSmartPointer<vtkVariantArray>::New());
299   crossEdge->InsertNextValue(1);
300   // CrossEdge is an edge linking hierarchies.
301   vtkUnsignedCharArray* crossEdgesArray=vtkUnsignedCharArray::New();
302   crossEdgesArray->SetName("CrossEdges");
303   sil->GetEdgeData()->AddArray(crossEdgesArray);
304   crossEdgesArray->Delete();
305   std::vector<std::string> names;
306   // Now build the hierarchy.
307   vtkIdType rootId=sil->AddVertex();
308   names.push_back("SIL");
309   // Add global fields root
310   vtkIdType fieldsRoot(sil->AddChild(rootId,childEdge));
311   names.push_back("FieldsStatusTree");
312   this->Internal->Tree.feedSIL(sil,fieldsRoot,childEdge,names);
313   vtkIdType meshesFamsGrpsRoot(sil->AddChild(rootId,childEdge));
314   names.push_back("MeshesFamsGrps");
315   std::string dftMeshName(this->Internal->Tree.feedSILForFamsAndGrps(sil,meshesFamsGrpsRoot,childEdge,names));
316   // This array is used to assign names to nodes.
317   vtkStringArray *namesArray(vtkStringArray::New());
318   namesArray->SetName("Names");
319   namesArray->SetNumberOfTuples(sil->GetNumberOfVertices());
320   sil->GetVertexData()->AddArray(namesArray);
321   namesArray->Delete();
322   std::vector<std::string>::const_iterator iter;
323   vtkIdType cc;
324   for(cc=0, iter=names.begin(); iter!=names.end(); ++iter, ++cc)
325     namesArray->SetValue(cc,(*iter).c_str());
326   return dftMeshName;
327 }
328
329 double vtkMEDReader::PublishTimeStepsIfNeeded(vtkInformation *outInfo, bool& isUpdated)
330 {
331   int lev0(-1);
332   std::vector<double> tsteps;
333   if(!this->Internal->IsStdOrMode)
334     tsteps=this->Internal->Tree.getTimeSteps(lev0,this->Internal->TK);
335   else
336     { tsteps.resize(1); tsteps[0]=0.; }
337   isUpdated=false;
338   if(lev0!=this->Internal->LastLev0)
339     {
340       isUpdated=true;
341       double timeRange[2];
342       timeRange[0]=tsteps.front();
343       timeRange[1]=tsteps.back();
344       outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),&tsteps[0],tsteps.size());
345       outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(),timeRange,2);
346       this->Internal->LastLev0=lev0;
347     }
348   return tsteps.front();
349 }
350
351 void vtkMEDReader::FillMultiBlockDataSetInstance(vtkMultiBlockDataSet *output, double reqTS)
352 {
353   std::string meshName;
354   vtkDataSet *ret(this->Internal->Tree.buildVTKInstance(this->Internal->IsStdOrMode,reqTS,meshName,this->Internal->TK));
355   if(this->Internal->GenerateVect)
356     {
357       vtkGenerateVectors::Operate(ret->GetPointData());
358       vtkGenerateVectors::Operate(ret->GetCellData());
359       vtkGenerateVectors::Operate(ret->GetFieldData());
360     }
361   output->SetBlock(0,ret);
362   ret->Delete();
363 }
364
365 void vtkMEDReader::PrintSelf(ostream& os, vtkIndent indent)
366 {
367   this->Superclass::PrintSelf(os, indent);
368 }