1 // Copyright (C) 2010-2022 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"
23 #include "MEDUtilities.hxx"
25 #include "vtkCellArray.h"
26 #include "vtkCellData.h"
27 #include "vtkCellType.h"
28 #include "vtkDataSetAttributes.h"
29 #include "vtkDataArraySelection.h"
30 #include "vtkDoubleArray.h"
31 #include "vtkExecutive.h"
32 #include "vtkInformation.h"
33 #include "vtkInformationDataObjectMetaDataKey.h"
34 #include "vtkInformationDoubleVectorKey.h"
35 #include "vtkInformationQuadratureSchemeDefinitionVectorKey.h"
36 #include "vtkInformationStringKey.h"
37 #include "vtkInformationVector.h"
38 #include "vtkMultiBlockDataSet.h"
39 #include "vtkMultiTimeStepAlgorithm.h"
40 #include "vtkMutableDirectedGraph.h"
41 #include "vtkObjectFactory.h"
42 #include "vtkPointData.h"
43 #include "vtkQuadratureSchemeDefinition.h"
44 #include "vtkSmartPointer.h"
45 #include "vtkStreamingDemandDrivenPipeline.h"
46 #include "vtkStringArray.h"
47 #include "vtkUnsignedCharArray.h"
48 #include "vtkUnstructuredGrid.h"
49 #include "vtkVariantArray.h"
51 #ifdef MEDREADER_USE_MPI
52 #include "vtkMultiProcessController.h"
53 #include "vtkGhostCellsGenerator.h"
56 #include "MEDFileFieldRepresentationTree.hxx"
64 class vtkMEDReader::vtkMEDReaderInternal
68 vtkMEDReaderInternal(vtkMEDReader *master):TK(0),SIL(0),LastLev0(-1)
72 ~vtkMEDReaderInternal()
78 MEDFileFieldRepresentationTree Tree;
81 std::string DftMeshName;
82 // Store the vtkMutableDirectedGraph that represents links between family, groups and cell types
83 vtkMutableDirectedGraph* SIL;
84 // store the lev0 id in Tree corresponding to the TIME_STEPS in the pipeline.
88 vtkStandardNewMacro(vtkMEDReader)
90 // vtkInformationKeyMacro(vtkMEDReader, META_DATA, DataObjectMetaData); // Here we need to customize vtkMEDReader::META_DATA method
91 // start of overload of vtkInformationKeyMacro
92 static vtkInformationDataObjectMetaDataKey *vtkMEDReader_META_DATA=new vtkInformationDataObjectMetaDataKey("META_DATA","vtkFileSeriesGroupReader");
94 vtkInformationDataObjectMetaDataKey *vtkMEDReader::META_DATA()
96 static const char ZE_KEY[]="vtkFileSeriesGroupReader::META_DATA";
97 vtkInformationDataObjectMetaDataKey *ret(vtkMEDReader_META_DATA);
98 MEDCoupling::GlobalDict *gd(MEDCoupling::GlobalDict::GetInstance());
99 if(!gd->hasKey(ZE_KEY))
100 {// here META_DATA is put on global var to be exchanged with other filters without dependancy of MEDReader. Please do not change ZE_KEY !
101 std::ostringstream oss; oss << ret;
102 gd->setKeyValue(ZE_KEY,oss.str());
107 static vtkInformationDoubleVectorKey *vtkMEDReader_GAUSS_DATA=new vtkInformationDoubleVectorKey("GAUSS_DATA","vtkFileSeriesGroupReader");
109 vtkInformationDoubleVectorKey *vtkMEDReader::GAUSS_DATA()
111 static const char ZE_KEY[]="vtkFileSeriesGroupReader::GAUSS_DATA";
112 vtkInformationDoubleVectorKey *ret(vtkMEDReader_GAUSS_DATA);
113 MEDCoupling::GlobalDict *gd(MEDCoupling::GlobalDict::GetInstance());
114 if(!gd->hasKey(ZE_KEY))
115 {// here META_DATA is put on global var to be exchanged with other filters without dependancy of MEDReader. Please do not change ZE_KEY !
116 vtkInformationDoubleVectorKey *ret2(ret);
117 std::ostringstream oss; oss << ret2;
118 gd->setKeyValue(ZE_KEY,oss.str());
122 // end of overload of vtkInformationKeyMacro
124 vtkMEDReader::vtkMEDReader():Internal(new vtkMEDReaderInternal(this))
126 this->SetNumberOfInputPorts(0);
127 this->SetNumberOfOutputPorts(1);
130 vtkMEDReader::~vtkMEDReader()
132 delete this->Internal;
136 void vtkMEDReader::Reload()
138 this->ReloadInternals();
139 this->IsStdOrMode = false;
140 this->GenerateVect = false;
142 this->FieldSelection->RemoveAllArrays();
143 this->TimeFlagSelection->RemoveAllArrays();
146 void vtkMEDReader::ReloadInternals()
148 delete this->Internal;
149 this->Internal=new vtkMEDReaderInternal(this);
153 void vtkMEDReader::GenerateVectors(int val)
155 if ( !this->Internal )
158 bool val2((bool)val);
159 if(val2!=this->GenerateVect)
161 this->GenerateVect=val2;
166 void vtkMEDReader::ChangeMode(int newMode)
168 if ( !this->Internal )
171 this->IsStdOrMode=newMode!=0;
175 void vtkMEDReader::GhostCellGeneratorCallForPara(int gcgcp)
177 if ( !this->Internal )
180 bool newVal(gcgcp!=0);
181 if(newVal!=this->GCGCP)
188 const char *vtkMEDReader::GetSeparator()
190 return MEDFileFieldRepresentationLeavesArrays::ZE_SEP;
193 void vtkMEDReader::SetFileName(const char *fname)
199 this->FileName=fname;
202 catch(INTERP_KERNEL::Exception& e)
204 delete this->Internal;
206 std::ostringstream oss;
207 oss << "Exception has been thrown in vtkMEDReader::SetFileName : " << e.what() << std::endl;
208 if(this->HasObserver("ErrorEvent") )
209 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
211 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
212 vtkObject::BreakOnError();
216 char *vtkMEDReader::GetFileName()
220 return const_cast<char *>(this->FileName.c_str());
223 int vtkMEDReader::RequestInformation(vtkInformation *request, vtkInformationVector ** /*inputVector*/, vtkInformationVector *outputVector)
225 // std::cout << "########################################## vtkMEDReader::RequestInformation ##########################################" << std::endl;
230 // Process file meta data
231 if(this->Internal->Tree.getNumberOfLeavesArrays()==0)
233 int iPart(-1),nbOfParts(-1);
234 #ifdef MEDREADER_USE_MPI
235 if (this->DistributeWithMPI)
237 vtkMultiProcessController *vmpc(vtkMultiProcessController::GetGlobalController());
240 iPart=vmpc->GetLocalProcessId();
241 nbOfParts=vmpc->GetNumberOfProcesses();
245 this->Internal->Tree.loadMainStructureOfFile(this->FileName.c_str(),iPart,nbOfParts);
248 this->Internal->Tree.activateTheFirst();//This line manually initialize the status of server (this) with the remote client.
249 for (int idLeaveArray = 0; idLeaveArray < this->Internal->Tree.getNumberOfLeavesArrays(); idLeaveArray++)
251 std::string name = this->Internal->Tree.getNameOf(idLeaveArray);
252 bool status = this->Internal->Tree.getStatusOf(idLeaveArray);
253 this->FieldSelection->AddArray(name.c_str(), status);
258 this->Internal->TK.setMaxNumberOfTimeSteps(this->Internal->Tree.getMaxNumberOfTimeSteps());
259 auto timeFlagsArray = this->Internal->TK.getTimesFlagArray();
260 for (int idTimeFlag = 0; idTimeFlag < timeFlagsArray.size() ; idTimeFlag++)
262 std::string name = timeFlagsArray[idTimeFlag].second;
263 bool status = timeFlagsArray[idTimeFlag].first;
264 this->TimeFlagSelection->AddArray(name.c_str(), status);
267 // Make sure internal model are synchronized
268 /// So the SIL is up to date
269 int nArrays = this->FieldSelection->GetNumberOfArrays();
270 for(int i = nArrays - 1; i >= 0; i--)
274 this->Internal->Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(
275 this->Internal->Tree.getIdHavingZeName(this->FieldSelection->GetArrayName(i)),
276 this->FieldSelection->GetArraySetting(i));
278 catch(INTERP_KERNEL::Exception& e)
280 // Remove the incorrect array
281 this->FieldSelection->RemoveArrayByIndex(i);
285 // request->Print(cout);
286 vtkInformation *outInfo(outputVector->GetInformationObject(0));
287 outInfo->Set(vtkDataObject::DATA_TYPE_NAME(),"vtkMultiBlockDataSet");
288 this->UpdateSIL(request, outInfo);
290 // Set the meta data graph as a meta data key in the information
291 // That's all that is needed to transfer it along the pipeline
292 outInfo->Set(vtkMEDReader::META_DATA(),this->Internal->SIL);
295 this->PublishTimeStepsIfNeeded(outInfo,dummy);
297 catch(INTERP_KERNEL::Exception& e)
299 std::ostringstream oss;
300 oss << "Exception has been thrown in vtkMEDReader::RequestInformation : " << e.what() << std::endl;
301 if(this->HasObserver("ErrorEvent") )
302 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
304 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
305 vtkObject::BreakOnError();
311 int vtkMEDReader::RequestData(vtkInformation *request, vtkInformationVector ** /*inputVector*/, vtkInformationVector *outputVector)
313 // std::cout << "########################################## vtkMEDReader::RequestData ##########################################" << std::endl;
318 for(int i = 0; i < this->FieldSelection->GetNumberOfArrays(); i++)
320 this->Internal->Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(
321 this->Internal->Tree.getIdHavingZeName(this->FieldSelection->GetArrayName(i)),
322 this->FieldSelection->GetArraySetting(i));
325 auto& timeFlagsArray = this->Internal->TK.getTimesFlagArray();
326 if (timeFlagsArray.size() != this->TimeFlagSelection->GetNumberOfArrays())
328 throw INTERP_KERNEL::Exception("Unexpected size of TimeFlagSelection");
330 for(int i = 0; i < this->TimeFlagSelection->GetNumberOfArrays(); i++)
332 timeFlagsArray[i] = std::make_pair(this->TimeFlagSelection->GetArraySetting(i),
333 this->TimeFlagSelection->GetArrayName(i));
336 // request->Print(cout);
337 vtkInformation *outInfo(outputVector->GetInformationObject(0));
338 vtkMultiBlockDataSet *output(vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
339 //bool isUpdated(false); // todo: unused
341 if(outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP()))
342 reqTS=outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP());
344 #ifndef MEDREADER_USE_MPI
345 this->FillMultiBlockDataSetInstance(output,reqTS,&ti);
347 if (this->DistributeWithMPI && this->GCGCP)
349 vtkSmartPointer<vtkGhostCellsGenerator> gcg(vtkSmartPointer<vtkGhostCellsGenerator>::New());
351 vtkDataSet *ret(RetrieveDataSetAtTime(reqTS,&ti));
352 gcg->SetInputData(ret);
356 // gcg->SetUseGlobalPointIds(true);
357 gcg->SetBuildIfRequired(false);
359 output->SetBlock(0,gcg->GetOutput());
362 this->FillMultiBlockDataSetInstance(output,reqTS,&ti);
366 const std::vector<double>& data(ti.getData());
367 outInfo->Set(vtkMEDReader::GAUSS_DATA(),&data[0],(int)data.size());
368 request->Append(vtkExecutive::KEYS_TO_COPY(),vtkMEDReader::GAUSS_DATA());// Thank you to SciberQuest and DIPOLE_CENTER ! Don't understand why ! In RequestInformation it does not work !
370 output->GetInformation()->Set(vtkDataObject::DATA_TIME_STEP(),reqTS);
371 // Is it really needed ? TODO
372 this->UpdateSIL(request, outInfo);
374 catch(INTERP_KERNEL::Exception& e)
376 std::cerr << "Exception has been thrown in vtkMEDReader::RequestData : " << e.what() << std::endl;
382 //------------------------------------------------------------------------------
383 int vtkMEDReader::GetNumberOfFieldsTreeArrays()
385 return this->FieldSelection->GetNumberOfArrays();
388 //------------------------------------------------------------------------------
389 const char* vtkMEDReader::GetFieldsTreeArrayName(int index)
391 return this->FieldSelection->GetArrayName(index);
394 //------------------------------------------------------------------------------
395 int vtkMEDReader::GetFieldsTreeArrayStatus(const char* name)
397 return this->FieldSelection->ArrayIsEnabled(name);
400 //------------------------------------------------------------------------------
401 void vtkMEDReader::SetFieldsStatus(const char* name, int status)
403 if (this->GetFieldsTreeArrayStatus(name) != status)
407 this->FieldSelection->EnableArray(name);
411 this->FieldSelection->DisableArray(name);
417 //------------------------------------------------------------------------------
418 int vtkMEDReader::GetNumberOfTimesFlagsArrays()
420 return this->TimeFlagSelection->GetNumberOfArrays();
423 //------------------------------------------------------------------------------
424 const char* vtkMEDReader::GetTimesFlagsArrayName(int index)
426 return this->TimeFlagSelection->GetArrayName(index);
429 //------------------------------------------------------------------------------
430 int vtkMEDReader::GetTimesFlagsArrayStatus(const char* name)
432 return this->TimeFlagSelection->ArrayIsEnabled(name);
435 //------------------------------------------------------------------------------
436 void vtkMEDReader::SetTimesFlagsStatus(const char* name, int status)
438 if (this->GetTimesFlagsArrayStatus(name) != status)
442 this->TimeFlagSelection->EnableArray(name);
446 this->TimeFlagSelection->DisableArray(name);
452 void vtkMEDReader::UpdateSIL(vtkInformation* /*request*/, vtkInformation * /*info*/)
456 std::string meshName(this->Internal->Tree.getActiveMeshName());
457 if(!this->Internal->SIL || meshName!=this->Internal->DftMeshName)
459 vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::New());
461 if(this->Internal->SIL)
462 this->Internal->SIL->Delete();
463 this->Internal->SIL=sil;
464 this->Internal->DftMeshName=meshName;
469 * The returned string is the name of the mesh activated which groups and families are in \a sil.
471 std::string vtkMEDReader::BuildSIL(vtkMutableDirectedGraph* sil)
474 return std::string();
476 vtkSmartPointer<vtkVariantArray> childEdge(vtkSmartPointer<vtkVariantArray>::New());
477 childEdge->InsertNextValue(0);
478 vtkSmartPointer<vtkVariantArray> crossEdge(vtkSmartPointer<vtkVariantArray>::New());
479 crossEdge->InsertNextValue(1);
480 // CrossEdge is an edge linking hierarchies.
481 vtkUnsignedCharArray* crossEdgesArray=vtkUnsignedCharArray::New();
482 crossEdgesArray->SetName("CrossEdges");
483 sil->GetEdgeData()->AddArray(crossEdgesArray);
484 crossEdgesArray->Delete();
485 std::vector<std::string> names;
486 // Now build the hierarchy.
487 vtkIdType rootId=sil->AddVertex();
488 names.push_back("SIL");
489 // Add global fields root
490 vtkIdType fieldsRoot(sil->AddChild(rootId,childEdge));
491 names.push_back("FieldsStatusTree");
492 this->Internal->Tree.feedSIL(sil,fieldsRoot,childEdge,names);
493 vtkIdType meshesFamsGrpsRoot(sil->AddChild(rootId,childEdge));
494 names.push_back("MeshesFamsGrps");
495 std::string dftMeshName(this->Internal->Tree.feedSILForFamsAndGrps(sil,meshesFamsGrpsRoot,childEdge,names));
496 // This array is used to assign names to nodes.
497 vtkStringArray *namesArray(vtkStringArray::New());
498 namesArray->SetName("Names");
499 namesArray->SetNumberOfTuples(sil->GetNumberOfVertices());
500 sil->GetVertexData()->AddArray(namesArray);
501 namesArray->Delete();
502 std::vector<std::string>::const_iterator iter;
504 for(cc=0, iter=names.begin(); iter!=names.end(); ++iter, ++cc)
505 namesArray->SetValue(cc,(*iter).c_str());
509 double vtkMEDReader::PublishTimeStepsIfNeeded(vtkInformation *outInfo, bool& isUpdated)
515 std::vector<double> tsteps;
516 if(!this->IsStdOrMode)
517 tsteps=this->Internal->Tree.getTimeSteps(lev0,this->Internal->TK);
519 { tsteps.resize(1); tsteps[0]=0.; }
521 if(lev0!=this->Internal->LastLev0)
525 timeRange[0]=tsteps.front();
526 timeRange[1]=tsteps.back();
527 outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),&tsteps[0],(int)tsteps.size());
528 outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(),timeRange,2);
529 this->Internal->LastLev0=lev0;
531 return tsteps.front();
534 void vtkMEDReader::FillMultiBlockDataSetInstance(vtkMultiBlockDataSet *output, double reqTS, ExportedTinyInfo *internalInfo)
536 if( !this->Internal )
538 vtkDataSet *ret(RetrieveDataSetAtTime(reqTS,internalInfo));
539 output->SetBlock(0,ret);
543 vtkDataSet *vtkMEDReader::RetrieveDataSetAtTime(double reqTS, ExportedTinyInfo *internalInfo)
545 if( !this->Internal )
547 std::string meshName;
548 vtkDataSet *ret(this->Internal->Tree.buildVTKInstance(this->IsStdOrMode,reqTS,meshName,this->Internal->TK,internalInfo));
549 if(this->GenerateVect)
551 vtkGenerateVectors::Operate(ret->GetPointData());
552 vtkGenerateVectors::Operate(ret->GetCellData());
553 vtkGenerateVectors::Operate(ret->GetFieldData());
554 // The operations above have potentially created new arrays -> This breaks the optimization of StaticMesh that expects the same field arrays over time.
555 // To enforce the cache recomputation declare modification of mesh.
556 //vtkGenerateVectors::ChangeMeshTimeToUpdateCache(ret);
561 void vtkMEDReader::PrintSelf(ostream& os, vtkIndent indent)
563 this->Superclass::PrintSelf(os, indent);