1 // Copyright (C) 2010-2016 CEA/DEN, EDF R&D
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.
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.
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
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
19 // Author : Anthony Geay
21 #include "vtkMEDReader.h"
22 #include "vtkGenerateVectors.h"
24 #include "vtkMultiBlockDataSet.h"
25 #include "vtkInformation.h"
26 #include "vtkDataSetAttributes.h"
27 #include "vtkStringArray.h"
28 #include "vtkMutableDirectedGraph.h"
29 #include "vtkInformationStringKey.h"
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 #include "vtkInformationDataObjectMetaDataKey.h"
49 #ifdef MEDREADER_USE_MPI
50 #include "vtkMultiProcessController.h"
53 #include "MEDFileFieldRepresentationTree.hxx"
62 * This class stores properties in loading state mode (pvsm) when the MED file has not been read yet.
63 * The file is not read beacause FileName has not been informed yet ! So this class stores properties of vtkMEDReader instance that
64 * owns it and wait the vtkMEDReader::SetFileName to apply properties afterwards.
69 PropertyKeeper(vtkMEDReader *master):_master(master),IsGVActivated(false),GVValue(0),IsCMActivated(false),CMValue(0) { }
70 void assignPropertiesIfNeeded();
71 bool arePropertiesOnTreeToSetAfter() const;
73 void pushFieldStatusEntry(const char* name, int status);
74 void pushGenerateVectorsValue(int value);
75 void pushChangeModeValue(int value);
76 void pushTimesFlagsStatusEntry(const char* name, int status);
78 // pool of pairs to assign in SetFieldsStatus if needed. The use case is the load using pvsm.
79 std::vector< std::pair<std::string,int> > SetFieldsStatusPairs;
87 std::vector< std::pair<std::string,int> > TimesFlagsStatusPairs;
88 vtkMEDReader *_master;
91 void PropertyKeeper::assignPropertiesIfNeeded()
93 if(!this->SetFieldsStatusPairs.empty())
95 for(std::vector< std::pair<std::string,int> >::const_iterator it=this->SetFieldsStatusPairs.begin();it!=this->SetFieldsStatusPairs.end();it++)
96 _master->SetFieldsStatus((*it).first.c_str(),(*it).second);
97 this->SetFieldsStatusPairs.clear();
99 if(!this->TimesFlagsStatusPairs.empty())
101 for(std::vector< std::pair<std::string,int> >::const_iterator it=this->TimesFlagsStatusPairs.begin();it!=this->TimesFlagsStatusPairs.end();it++)
102 _master->SetTimesFlagsStatus((*it).first.c_str(),(*it).second);
103 this->TimesFlagsStatusPairs.clear();
105 if(this->IsGVActivated)
107 _master->GenerateVectors(this->GVValue);
108 this->IsGVActivated=false;
110 if(this->IsCMActivated)
112 _master->ChangeMode(this->CMValue);
113 this->IsCMActivated=false;
117 void PropertyKeeper::pushFieldStatusEntry(const char* name, int status)
120 for(std::vector< std::pair<std::string,int> >::const_iterator it=this->SetFieldsStatusPairs.begin();it!=this->SetFieldsStatusPairs.end() && !found;it++)
121 found=(*it).first==name;
123 this->SetFieldsStatusPairs.push_back(std::pair<std::string,int>(name,status));
126 void PropertyKeeper::pushTimesFlagsStatusEntry(const char* name, int status)
129 for(std::vector< std::pair<std::string,int> >::const_iterator it=this->TimesFlagsStatusPairs.begin();it!=this->TimesFlagsStatusPairs.end() && !found;it++)
130 found=(*it).first==name;
132 this->TimesFlagsStatusPairs.push_back(std::pair<std::string,int>(name,status));
135 void PropertyKeeper::pushGenerateVectorsValue(int value)
137 this->IsGVActivated=true;
141 void PropertyKeeper::pushChangeModeValue(int value)
143 this->IsCMActivated=true;
147 bool PropertyKeeper::arePropertiesOnTreeToSetAfter() const
149 return !SetFieldsStatusPairs.empty();
152 class vtkMEDReader::vtkMEDReaderInternal
156 vtkMEDReaderInternal(vtkMEDReader *master):TK(0),IsMEDOrSauv(true),IsStdOrMode(false),GenerateVect(false),SIL(0),LastLev0(-1),FirstCall0(2),PK(master),MyMTime(0)
162 return false; // TODO Useless and buggy
169 ~vtkMEDReaderInternal()
175 MEDFileFieldRepresentationTree Tree;
177 std::string FileName;
178 //when true the file is MED file. when false it is a Sauv file
180 //when false -> std, true -> mode. By default std (false).
182 //when false -> do nothing. When true cut off or extend to nbOfCompo=3 vector arrays.
184 std::string DftMeshName;
185 // Store the vtkMutableDirectedGraph that represents links between family, groups and cell types
186 vtkMutableDirectedGraph* SIL;
187 // store the lev0 id in Tree corresponding to the TIME_STEPS in the pipeline.
189 // The property keeper is usable only in pvsm mode.
192 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 ?
193 std::map<std::string,bool> _wonderful_ref;// this map stores the state before a SetFieldsStatus status.
196 unsigned char FirstCall0;
199 vtkStandardNewMacro(vtkMEDReader);
201 // vtkInformationKeyMacro(vtkMEDReader, META_DATA, DataObjectMetaData); // Here we need to customize vtkMEDReader::META_DATA method
202 // start of overload of vtkInformationKeyMacro
203 static vtkInformationDataObjectMetaDataKey *vtkMEDReader_META_DATA=new vtkInformationDataObjectMetaDataKey("META_DATA","vtkMEDReader");
205 vtkInformationDataObjectMetaDataKey *vtkMEDReader::META_DATA()
207 static const char ZE_KEY[]="vtkMEDReader::META_DATA";
208 vtkInformationDataObjectMetaDataKey *ret(vtkMEDReader_META_DATA);
209 MEDCoupling::GlobalDict *gd(MEDCoupling::GlobalDict::GetInstance());
210 if(!gd->hasKey(ZE_KEY))
211 {// here META_DATA is put on global var to be exchanged with other filters without dependancy of MEDReader. Please do not change ZE_KEY !
212 std::ostringstream oss; oss << ret;
213 gd->setKeyValue(ZE_KEY,oss.str());
217 // end of overload of vtkInformationKeyMacro
219 vtkMEDReader::vtkMEDReader():Internal(new vtkMEDReaderInternal(this))
221 this->SetNumberOfInputPorts(0);
222 this->SetNumberOfOutputPorts(1);
225 vtkMEDReader::~vtkMEDReader()
227 delete this->Internal;
231 void vtkMEDReader::Reload()
233 std::string fName((const char *)this->GetFileName());
234 delete this->Internal;
235 this->Internal=new vtkMEDReaderInternal(this);
236 this->SetFileName(fName.c_str());
239 int vtkMEDReader::GetServerModifTime()
241 if( !this->Internal )
243 return this->Internal->MyMTime;
246 void vtkMEDReader::GenerateVectors(int val)
248 if ( !this->Internal )
251 if(this->Internal->FileName.empty())
253 this->Internal->PK.pushGenerateVectorsValue(val);
256 //not pvsm mode (general case)
257 bool val2((bool)val);
258 if(val2!=this->Internal->GenerateVect)
260 this->Internal->GenerateVect=val2;
265 void vtkMEDReader::ChangeMode(int newMode)
267 if ( !this->Internal )
270 if(this->Internal->FileName.empty())
272 this->Internal->PK.pushChangeModeValue(newMode);
275 //not pvsm mode (general case)
276 this->Internal->IsStdOrMode=newMode!=0;
280 const char *vtkMEDReader::GetSeparator()
282 return MEDFileFieldRepresentationLeavesArrays::ZE_SEP;
285 void vtkMEDReader::SetFileName(const char *fname)
291 this->Internal->FileName=fname;
292 std::size_t pos(this->Internal->FileName.find_last_of('.'));
293 if(pos!=std::string::npos)
295 std::string ext(this->Internal->FileName.substr(pos));
296 if(ext.find("sauv")!=std::string::npos)
297 this->Internal->IsMEDOrSauv=false;
299 if(this->Internal->Tree.getNumberOfLeavesArrays()==0)
301 int iPart(-1),nbOfParts(-1);
302 #ifdef MEDREADER_USE_MPI
303 vtkMultiProcessController *vmpc(vtkMultiProcessController::GetGlobalController());
306 iPart=vmpc->GetLocalProcessId();
307 nbOfParts=vmpc->GetNumberOfProcesses();
310 this->Internal->Tree.loadMainStructureOfFile(this->Internal->FileName.c_str(),this->Internal->IsMEDOrSauv,iPart,nbOfParts);
311 if(!this->Internal->PK.arePropertiesOnTreeToSetAfter())
312 this->Internal->Tree.activateTheFirst();//This line manually initialize the status of server (this) with the remote client.
313 this->Internal->TK.setMaxNumberOfTimeSteps(this->Internal->Tree.getMaxNumberOfTimeSteps());
316 this->Internal->PK.assignPropertiesIfNeeded();
318 catch(INTERP_KERNEL::Exception& e)
320 delete this->Internal;
322 std::ostringstream oss;
323 oss << "Exception has been thrown in vtkMEDReader::SetFileName : " << e.what() << std::endl;
324 if(this->HasObserver("ErrorEvent") )
325 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
327 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
328 vtkObject::BreakOnError();
332 char *vtkMEDReader::GetFileName()
336 return const_cast<char *>(this->Internal->FileName.c_str());
339 int vtkMEDReader::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
341 // std::cout << "########################################## vtkMEDReader::RequestInformation ##########################################" << std::endl;
346 // request->Print(cout);
347 vtkInformation *outInfo(outputVector->GetInformationObject(0));
348 outInfo->Set(vtkDataObject::DATA_TYPE_NAME(),"vtkMultiBlockDataSet");
349 this->UpdateSIL(request, outInfo);
351 // Set the meta data graph as a meta data key in the information
352 // That's all that is needed to transfer it along the pipeline
353 outInfo->Set(vtkMEDReader::META_DATA(),this->Internal->SIL);
356 this->PublishTimeStepsIfNeeded(outInfo,dummy);
358 catch(INTERP_KERNEL::Exception& e)
360 std::ostringstream oss;
361 oss << "Exception has been thrown in vtkMEDReader::RequestInformation : " << e.what() << std::endl;
362 if(this->HasObserver("ErrorEvent") )
363 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
365 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
366 vtkObject::BreakOnError();
372 int vtkMEDReader::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
374 // std::cout << "########################################## vtkMEDReader::RequestData ##########################################" << std::endl;
379 // request->Print(cout);
380 vtkInformation *outInfo(outputVector->GetInformationObject(0));
381 vtkMultiBlockDataSet *output(vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
382 bool isUpdated(false);
384 if(outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP()))
385 reqTS=outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP());
386 this->FillMultiBlockDataSetInstance(output,reqTS);
387 output->GetInformation()->Set(vtkDataObject::DATA_TIME_STEP(),reqTS);
389 // Is it really needed ? TODO
390 this->UpdateSIL(request, outInfo);
392 catch(INTERP_KERNEL::Exception& e)
394 std::cerr << "Exception has been thrown in vtkMEDReader::RequestData : " << e.what() << std::endl;
400 void vtkMEDReader::SetFieldsStatus(const char* name, int status)
402 if( !this->Internal )
405 //this->Internal->_wonderful_set.insert(name);
406 if(this->Internal->FileName.empty())
408 this->Internal->PK.pushFieldStatusEntry(name,status);
411 if(this->Internal->_wonderful_set.empty())
412 this->Internal->_wonderful_ref=this->Internal->Tree.dumpState();// start of SetFieldsStatus serie -> store ref to compare at the end of the SetFieldsStatus serie.
413 this->Internal->_wonderful_set.insert(name);
414 //not pvsm mode (general case)
417 this->Internal->Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(this->Internal->Tree.getIdHavingZeName(name),status);
418 if(this->Internal->_wonderful_set.size()==GetNumberOfFieldsTreeArrays())
420 if(this->Internal->_wonderful_ref!=this->Internal->Tree.dumpState())
422 if(!this->Internal->PluginStart0())
426 this->Internal->MyMTime++;
428 this->Internal->_wonderful_set.clear();
431 catch(INTERP_KERNEL::Exception& e)
433 if(!this->Internal->FileName.empty())
435 std::cerr << "vtkMEDReader::SetFieldsStatus error : " << e.what() << " *** WITH STATUS=" << status << endl;
441 int vtkMEDReader::GetNumberOfFieldsTreeArrays()
445 return this->Internal->Tree.getNumberOfLeavesArrays();
448 const char *vtkMEDReader::GetFieldsTreeArrayName(int index)
452 return this->Internal->Tree.getNameOfC(index);
455 int vtkMEDReader::GetFieldsTreeArrayStatus(const char *name)
460 int zeId(this->Internal->Tree.getIdHavingZeName(name));
461 int ret(this->Internal->Tree.getStatusOf(zeId));
465 void vtkMEDReader::SetTimesFlagsStatus(const char *name, int status)
470 if(this->Internal->FileName.empty())
472 this->Internal->PK.pushTimesFlagsStatusEntry(name,status);
475 //not pvsm mode (general case)
477 std::istringstream iss(name); iss >> pos;
478 this->Internal->TK.getTimesFlagArray()[pos].first=(bool)status;
479 if(pos==this->Internal->TK.getTimesFlagArray().size()-1)
480 if(!this->Internal->PluginStart0())
483 //this->Internal->TK.printSelf(std::cerr);
487 int vtkMEDReader::GetNumberOfTimesFlagsArrays()
491 return (int)this->Internal->TK.getTimesFlagArray().size();
494 const char *vtkMEDReader::GetTimesFlagsArrayName(int index)
496 return this->Internal->TK.getTimesFlagArray()[index].second.c_str();
499 int vtkMEDReader::GetTimesFlagsArrayStatus(const char *name)
504 std::istringstream iss(name); iss >> pos;
505 return (int)this->Internal->TK.getTimesFlagArray()[pos].first;
508 void vtkMEDReader::UpdateSIL(vtkInformation* request, vtkInformation *info)
512 std::string meshName(this->Internal->Tree.getActiveMeshName());
513 if(!this->Internal->SIL || meshName!=this->Internal->DftMeshName)
515 vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::New());
517 if(this->Internal->SIL)
518 this->Internal->SIL->Delete();
519 this->Internal->SIL=sil;
520 this->Internal->DftMeshName=meshName;
525 * The returned string is the name of the mesh activated which groups and families are in \a sil.
527 std::string vtkMEDReader::BuildSIL(vtkMutableDirectedGraph* sil)
530 return std::string();
532 vtkSmartPointer<vtkVariantArray> childEdge(vtkSmartPointer<vtkVariantArray>::New());
533 childEdge->InsertNextValue(0);
534 vtkSmartPointer<vtkVariantArray> crossEdge(vtkSmartPointer<vtkVariantArray>::New());
535 crossEdge->InsertNextValue(1);
536 // CrossEdge is an edge linking hierarchies.
537 vtkUnsignedCharArray* crossEdgesArray=vtkUnsignedCharArray::New();
538 crossEdgesArray->SetName("CrossEdges");
539 sil->GetEdgeData()->AddArray(crossEdgesArray);
540 crossEdgesArray->Delete();
541 std::vector<std::string> names;
542 // Now build the hierarchy.
543 vtkIdType rootId=sil->AddVertex();
544 names.push_back("SIL");
545 // Add global fields root
546 vtkIdType fieldsRoot(sil->AddChild(rootId,childEdge));
547 names.push_back("FieldsStatusTree");
548 this->Internal->Tree.feedSIL(sil,fieldsRoot,childEdge,names);
549 vtkIdType meshesFamsGrpsRoot(sil->AddChild(rootId,childEdge));
550 names.push_back("MeshesFamsGrps");
551 std::string dftMeshName(this->Internal->Tree.feedSILForFamsAndGrps(sil,meshesFamsGrpsRoot,childEdge,names));
552 // This array is used to assign names to nodes.
553 vtkStringArray *namesArray(vtkStringArray::New());
554 namesArray->SetName("Names");
555 namesArray->SetNumberOfTuples(sil->GetNumberOfVertices());
556 sil->GetVertexData()->AddArray(namesArray);
557 namesArray->Delete();
558 std::vector<std::string>::const_iterator iter;
560 for(cc=0, iter=names.begin(); iter!=names.end(); ++iter, ++cc)
561 namesArray->SetValue(cc,(*iter).c_str());
565 double vtkMEDReader::PublishTimeStepsIfNeeded(vtkInformation *outInfo, bool& isUpdated)
571 std::vector<double> tsteps;
572 if(!this->Internal->IsStdOrMode)
573 tsteps=this->Internal->Tree.getTimeSteps(lev0,this->Internal->TK);
575 { tsteps.resize(1); tsteps[0]=0.; }
577 if(lev0!=this->Internal->LastLev0)
581 timeRange[0]=tsteps.front();
582 timeRange[1]=tsteps.back();
583 outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),&tsteps[0],tsteps.size());
584 outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(),timeRange,2);
585 this->Internal->LastLev0=lev0;
587 return tsteps.front();
590 void vtkMEDReader::FillMultiBlockDataSetInstance(vtkMultiBlockDataSet *output, double reqTS)
592 if( !this->Internal )
594 std::string meshName;
595 vtkDataSet *ret(this->Internal->Tree.buildVTKInstance(this->Internal->IsStdOrMode,reqTS,meshName,this->Internal->TK));
596 if(this->Internal->GenerateVect)
598 vtkGenerateVectors::Operate(ret->GetPointData());
599 vtkGenerateVectors::Operate(ret->GetCellData());
600 vtkGenerateVectors::Operate(ret->GetFieldData());
602 output->SetBlock(0,ret);
606 void vtkMEDReader::PrintSelf(ostream& os, vtkIndent indent)
608 this->Superclass::PrintSelf(os, indent);