1 // Copyright (C) 2010-2021 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 "vtkPUnstructuredGridGhostCellsGenerator.h"
56 #include "MEDFileFieldRepresentationTree.hxx"
64 class vtkMEDReader::vtkMEDReaderInternal
68 vtkMEDReaderInternal(vtkMEDReader *master):TK(0),IsStdOrMode(false),GenerateVect(false),SIL(0),LastLev0(-1),GCGCP(true)
72 ~vtkMEDReaderInternal()
78 MEDFileFieldRepresentationTree Tree;
79 vtkNew<vtkDataArraySelection> FieldSelection;
80 vtkNew<vtkDataArraySelection> TimeFlagSelection;
84 //when false -> std, true -> mode. By default std (false).
86 //when false -> do nothing. When true cut off or extend to nbOfCompo=3 vector arrays.
88 std::string DftMeshName;
89 // Store the vtkMutableDirectedGraph that represents links between family, groups and cell types
90 vtkMutableDirectedGraph* SIL;
91 // store the lev0 id in Tree corresponding to the TIME_STEPS in the pipeline.
96 vtkStandardNewMacro(vtkMEDReader)
98 // vtkInformationKeyMacro(vtkMEDReader, META_DATA, DataObjectMetaData); // Here we need to customize vtkMEDReader::META_DATA method
99 // start of overload of vtkInformationKeyMacro
100 static vtkInformationDataObjectMetaDataKey *vtkMEDReader_META_DATA=new vtkInformationDataObjectMetaDataKey("META_DATA","vtkMEDReader");
102 vtkInformationDataObjectMetaDataKey *vtkMEDReader::META_DATA()
104 static const char ZE_KEY[]="vtkMEDReader::META_DATA";
105 vtkInformationDataObjectMetaDataKey *ret(vtkMEDReader_META_DATA);
106 MEDCoupling::GlobalDict *gd(MEDCoupling::GlobalDict::GetInstance());
107 if(!gd->hasKey(ZE_KEY))
108 {// here META_DATA is put on global var to be exchanged with other filters without dependancy of MEDReader. Please do not change ZE_KEY !
109 std::ostringstream oss; oss << ret;
110 gd->setKeyValue(ZE_KEY,oss.str());
115 static vtkInformationGaussDoubleVectorKey *vtkMEDReader_GAUSS_DATA=new vtkInformationGaussDoubleVectorKey("GAUSS_DATA","vtkMEDReader");
117 vtkInformationGaussDoubleVectorKey *vtkMEDReader::GAUSS_DATA()
119 static const char ZE_KEY[]="vtkMEDReader::GAUSS_DATA";
120 vtkInformationGaussDoubleVectorKey *ret(vtkMEDReader_GAUSS_DATA);
121 MEDCoupling::GlobalDict *gd(MEDCoupling::GlobalDict::GetInstance());
122 if(!gd->hasKey(ZE_KEY))
123 {// here META_DATA is put on global var to be exchanged with other filters without dependancy of MEDReader. Please do not change ZE_KEY !
124 vtkInformationDoubleVectorKey *ret2(ret);
125 std::ostringstream oss; oss << ret2;
126 gd->setKeyValue(ZE_KEY,oss.str());
130 // end of overload of vtkInformationKeyMacro
132 vtkMEDReader::vtkMEDReader():Internal(new vtkMEDReaderInternal(this))
134 this->SetNumberOfInputPorts(0);
135 this->SetNumberOfOutputPorts(1);
138 vtkMEDReader::~vtkMEDReader()
140 delete this->Internal;
144 void vtkMEDReader::Reload()
146 std::string fName((const char *)this->GetFileName());
147 delete this->Internal;
148 this->Internal=new vtkMEDReaderInternal(this);
149 this->SetFileName(fName.c_str());
152 void vtkMEDReader::GenerateVectors(int val)
154 if ( !this->Internal )
157 bool val2((bool)val);
158 if(val2!=this->Internal->GenerateVect)
160 this->Internal->GenerateVect=val2;
165 void vtkMEDReader::ChangeMode(int newMode)
167 if ( !this->Internal )
170 this->Internal->IsStdOrMode=newMode!=0;
174 void vtkMEDReader::GhostCellGeneratorCallForPara(int gcgcp)
176 if ( !this->Internal )
179 bool newVal(gcgcp!=0);
180 if(newVal!=this->Internal->GCGCP)
182 this->Internal->GCGCP=newVal;
187 const char *vtkMEDReader::GetSeparator()
189 return MEDFileFieldRepresentationLeavesArrays::ZE_SEP;
192 void vtkMEDReader::SetFileName(const char *fname)
198 this->Internal->FileName=fname;
201 catch(INTERP_KERNEL::Exception& e)
203 delete this->Internal;
205 std::ostringstream oss;
206 oss << "Exception has been thrown in vtkMEDReader::SetFileName : " << e.what() << std::endl;
207 if(this->HasObserver("ErrorEvent") )
208 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
210 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
211 vtkObject::BreakOnError();
215 char *vtkMEDReader::GetFileName()
219 return const_cast<char *>(this->Internal->FileName.c_str());
222 int vtkMEDReader::RequestInformation(vtkInformation *request, vtkInformationVector ** /*inputVector*/, vtkInformationVector *outputVector)
224 // std::cout << "########################################## vtkMEDReader::RequestInformation ##########################################" << std::endl;
229 // Process file meta data
230 if(this->Internal->Tree.getNumberOfLeavesArrays()==0)
232 int iPart(-1),nbOfParts(-1);
233 #ifdef MEDREADER_USE_MPI
234 vtkMultiProcessController *vmpc(vtkMultiProcessController::GetGlobalController());
237 iPart=vmpc->GetLocalProcessId();
238 nbOfParts=vmpc->GetNumberOfProcesses();
241 this->Internal->Tree.loadMainStructureOfFile(this->Internal->FileName.c_str(),iPart,nbOfParts);
244 this->Internal->Tree.activateTheFirst();//This line manually initialize the status of server (this) with the remote client.
245 for (int idLeaveArray = 0; idLeaveArray < this->Internal->Tree.getNumberOfLeavesArrays(); idLeaveArray++)
247 std::string name = this->Internal->Tree.getNameOf(idLeaveArray);
248 bool status = this->Internal->Tree.getStatusOf(idLeaveArray);
249 this->Internal->FieldSelection->AddArray(name.c_str(), status);
254 this->Internal->TK.setMaxNumberOfTimeSteps(this->Internal->Tree.getMaxNumberOfTimeSteps());
255 auto timeFlagsArray = this->Internal->TK.getTimesFlagArray();
256 for (int idTimeFlag = 0; idTimeFlag < timeFlagsArray.size() ; idTimeFlag++)
258 std::string name = timeFlagsArray[idTimeFlag].second;
259 bool status = timeFlagsArray[idTimeFlag].first;
260 this->Internal->TimeFlagSelection->AddArray(name.c_str(), status);
263 // Make sure internal model are synchronized
264 /// So the SIL is up to date
265 int nArrays = this->Internal->FieldSelection->GetNumberOfArrays();
266 for(int i = nArrays - 1; i >= 0; i--)
270 this->Internal->Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(
271 this->Internal->Tree.getIdHavingZeName(this->Internal->FieldSelection->GetArrayName(i)),
272 this->Internal->FieldSelection->GetArraySetting(i));
274 catch(INTERP_KERNEL::Exception& e)
276 // Remove the incorrect array
277 this->Internal->FieldSelection->RemoveArrayByIndex(i);
281 // request->Print(cout);
282 vtkInformation *outInfo(outputVector->GetInformationObject(0));
283 outInfo->Set(vtkDataObject::DATA_TYPE_NAME(),"vtkMultiBlockDataSet");
284 this->UpdateSIL(request, outInfo);
286 // Set the meta data graph as a meta data key in the information
287 // That's all that is needed to transfer it along the pipeline
288 outInfo->Set(vtkMEDReader::META_DATA(),this->Internal->SIL);
291 this->PublishTimeStepsIfNeeded(outInfo,dummy);
293 catch(INTERP_KERNEL::Exception& e)
295 std::ostringstream oss;
296 oss << "Exception has been thrown in vtkMEDReader::RequestInformation : " << e.what() << std::endl;
297 if(this->HasObserver("ErrorEvent") )
298 this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
300 vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
301 vtkObject::BreakOnError();
307 int vtkMEDReader::RequestData(vtkInformation *request, vtkInformationVector ** /*inputVector*/, vtkInformationVector *outputVector)
309 // std::cout << "########################################## vtkMEDReader::RequestData ##########################################" << std::endl;
314 for(int i = 0; i < this->Internal->FieldSelection->GetNumberOfArrays(); i++)
316 this->Internal->Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(
317 this->Internal->Tree.getIdHavingZeName(this->Internal->FieldSelection->GetArrayName(i)),
318 this->Internal->FieldSelection->GetArraySetting(i));
321 auto& timeFlagsArray = this->Internal->TK.getTimesFlagArray();
322 if (timeFlagsArray.size() != this->Internal->TimeFlagSelection->GetNumberOfArrays())
324 throw INTERP_KERNEL::Exception("Unexpected size of TimeFlagSelection");
326 for(int i = 0; i < this->Internal->TimeFlagSelection->GetNumberOfArrays(); i++)
328 timeFlagsArray[i] = std::make_pair(this->Internal->TimeFlagSelection->GetArraySetting(i),
329 this->Internal->TimeFlagSelection->GetArrayName(i));
332 // request->Print(cout);
333 vtkInformation *outInfo(outputVector->GetInformationObject(0));
334 vtkMultiBlockDataSet *output(vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
335 //bool isUpdated(false); // todo: unused
337 if(outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP()))
338 reqTS=outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP());
340 #ifndef MEDREADER_USE_MPI
341 this->FillMultiBlockDataSetInstance(output,reqTS,&ti);
343 if(this->Internal->GCGCP)
345 vtkSmartPointer<vtkPUnstructuredGridGhostCellsGenerator> gcg(vtkSmartPointer<vtkPUnstructuredGridGhostCellsGenerator>::New());
347 vtkDataSet *ret(RetrieveDataSetAtTime(reqTS,&ti));
348 gcg->SetInputData(ret);
351 gcg->SetUseGlobalPointIds(true);
352 gcg->SetBuildIfRequired(false);
354 output->SetBlock(0,gcg->GetOutput());
357 this->FillMultiBlockDataSetInstance(output,reqTS,&ti);
361 const std::vector<double>& data(ti.getData());
362 outInfo->Set(vtkMEDReader::GAUSS_DATA(),&data[0],(int)data.size());
363 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 !
365 output->GetInformation()->Set(vtkDataObject::DATA_TIME_STEP(),reqTS);
366 // Is it really needed ? TODO
367 this->UpdateSIL(request, outInfo);
369 catch(INTERP_KERNEL::Exception& e)
371 std::cerr << "Exception has been thrown in vtkMEDReader::RequestData : " << e.what() << std::endl;
377 //------------------------------------------------------------------------------
378 int vtkMEDReader::GetNumberOfFieldsTreeArrays()
380 return this->Internal->FieldSelection->GetNumberOfArrays();
383 //------------------------------------------------------------------------------
384 const char* vtkMEDReader::GetFieldsTreeArrayName(int index)
386 return this->Internal->FieldSelection->GetArrayName(index);
389 //------------------------------------------------------------------------------
390 int vtkMEDReader::GetFieldsTreeArrayStatus(const char* name)
392 return this->Internal->FieldSelection->ArrayIsEnabled(name);
395 //------------------------------------------------------------------------------
396 void vtkMEDReader::SetFieldsStatus(const char* name, int status)
398 if (this->GetFieldsTreeArrayStatus(name) != status)
402 this->Internal->FieldSelection->EnableArray(name);
406 this->Internal->FieldSelection->DisableArray(name);
412 //------------------------------------------------------------------------------
413 int vtkMEDReader::GetNumberOfTimesFlagsArrays()
415 return this->Internal->TimeFlagSelection->GetNumberOfArrays();
418 //------------------------------------------------------------------------------
419 const char* vtkMEDReader::GetTimesFlagsArrayName(int index)
421 return this->Internal->TimeFlagSelection->GetArrayName(index);
424 //------------------------------------------------------------------------------
425 int vtkMEDReader::GetTimesFlagsArrayStatus(const char* name)
427 return this->Internal->TimeFlagSelection->ArrayIsEnabled(name);
430 //------------------------------------------------------------------------------
431 void vtkMEDReader::SetTimesFlagsStatus(const char* name, int status)
433 if (this->GetTimesFlagsArrayStatus(name) != status)
437 this->Internal->TimeFlagSelection->EnableArray(name);
441 this->Internal->TimeFlagSelection->DisableArray(name);
447 void vtkMEDReader::UpdateSIL(vtkInformation* /*request*/, vtkInformation * /*info*/)
451 std::string meshName(this->Internal->Tree.getActiveMeshName());
452 if(!this->Internal->SIL || meshName!=this->Internal->DftMeshName)
454 vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::New());
456 if(this->Internal->SIL)
457 this->Internal->SIL->Delete();
458 this->Internal->SIL=sil;
459 this->Internal->DftMeshName=meshName;
464 * The returned string is the name of the mesh activated which groups and families are in \a sil.
466 std::string vtkMEDReader::BuildSIL(vtkMutableDirectedGraph* sil)
469 return std::string();
471 vtkSmartPointer<vtkVariantArray> childEdge(vtkSmartPointer<vtkVariantArray>::New());
472 childEdge->InsertNextValue(0);
473 vtkSmartPointer<vtkVariantArray> crossEdge(vtkSmartPointer<vtkVariantArray>::New());
474 crossEdge->InsertNextValue(1);
475 // CrossEdge is an edge linking hierarchies.
476 vtkUnsignedCharArray* crossEdgesArray=vtkUnsignedCharArray::New();
477 crossEdgesArray->SetName("CrossEdges");
478 sil->GetEdgeData()->AddArray(crossEdgesArray);
479 crossEdgesArray->Delete();
480 std::vector<std::string> names;
481 // Now build the hierarchy.
482 vtkIdType rootId=sil->AddVertex();
483 names.push_back("SIL");
484 // Add global fields root
485 vtkIdType fieldsRoot(sil->AddChild(rootId,childEdge));
486 names.push_back("FieldsStatusTree");
487 this->Internal->Tree.feedSIL(sil,fieldsRoot,childEdge,names);
488 vtkIdType meshesFamsGrpsRoot(sil->AddChild(rootId,childEdge));
489 names.push_back("MeshesFamsGrps");
490 std::string dftMeshName(this->Internal->Tree.feedSILForFamsAndGrps(sil,meshesFamsGrpsRoot,childEdge,names));
491 // This array is used to assign names to nodes.
492 vtkStringArray *namesArray(vtkStringArray::New());
493 namesArray->SetName("Names");
494 namesArray->SetNumberOfTuples(sil->GetNumberOfVertices());
495 sil->GetVertexData()->AddArray(namesArray);
496 namesArray->Delete();
497 std::vector<std::string>::const_iterator iter;
499 for(cc=0, iter=names.begin(); iter!=names.end(); ++iter, ++cc)
500 namesArray->SetValue(cc,(*iter).c_str());
504 double vtkMEDReader::PublishTimeStepsIfNeeded(vtkInformation *outInfo, bool& isUpdated)
510 std::vector<double> tsteps;
511 if(!this->Internal->IsStdOrMode)
512 tsteps=this->Internal->Tree.getTimeSteps(lev0,this->Internal->TK);
514 { tsteps.resize(1); tsteps[0]=0.; }
516 if(lev0!=this->Internal->LastLev0)
520 timeRange[0]=tsteps.front();
521 timeRange[1]=tsteps.back();
522 outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),&tsteps[0],(int)tsteps.size());
523 outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(),timeRange,2);
524 this->Internal->LastLev0=lev0;
526 return tsteps.front();
529 void vtkMEDReader::FillMultiBlockDataSetInstance(vtkMultiBlockDataSet *output, double reqTS, ExportedTinyInfo *internalInfo)
531 if( !this->Internal )
533 vtkDataSet *ret(RetrieveDataSetAtTime(reqTS,internalInfo));
534 output->SetBlock(0,ret);
538 vtkDataSet *vtkMEDReader::RetrieveDataSetAtTime(double reqTS, ExportedTinyInfo *internalInfo)
540 if( !this->Internal )
542 std::string meshName;
543 vtkDataSet *ret(this->Internal->Tree.buildVTKInstance(this->Internal->IsStdOrMode,reqTS,meshName,this->Internal->TK,internalInfo));
544 if(this->Internal->GenerateVect)
546 vtkGenerateVectors::Operate(ret->GetPointData());
547 vtkGenerateVectors::Operate(ret->GetCellData());
548 vtkGenerateVectors::Operate(ret->GetFieldData());
549 // The operations above have potentially created new arrays -> This breaks the optimization of StaticMesh that expects the same field arrays over time.
550 // To enforce the cache recomputation declare modification of mesh.
551 //vtkGenerateVectors::ChangeMeshTimeToUpdateCache(ret);
556 void vtkMEDReader::PrintSelf(ostream& os, vtkIndent indent)
558 this->Superclass::PrintSelf(os, indent);