From: Ovidiu Mircescu Date: Thu, 11 Aug 2016 09:27:51 +0000 (+0200) Subject: Add Catalyst processing to Insitu library. X-Git-Tag: SHAPER_2.7.0~15^2~6 X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=3c3ee02146d131a53e0567c32b0e111cf66fa6b3;p=modules%2Fparavis.git Add Catalyst processing to Insitu library. --- diff --git a/src/Insitu/CMakeLists.txt b/src/Insitu/CMakeLists.txt index 63c1169e..62d1db06 100644 --- a/src/Insitu/CMakeLists.txt +++ b/src/Insitu/CMakeLists.txt @@ -37,6 +37,6 @@ INSTALL(TARGETS visulib LIBRARY DESTINATION lib/salome ARCHIVE DESTINATION lib/salome ) -INSTALL(FILES ${_lib_HEADERS} DESTINATION include) +INSTALL(FILES ${_lib_HEADERS} DESTINATION include/salome ) diff --git a/src/Insitu/visu.cxx b/src/Insitu/visu.cxx old mode 100644 new mode 100755 index a454d80d..9c85828f --- a/src/Insitu/visu.cxx +++ b/src/Insitu/visu.cxx @@ -1,5 +1,4 @@ #include "visu.hxx" -// #include "MEDCouplingFieldDouble.hxx" #include "MEDFileField.hxx" #include "MEDCouplingUMesh.hxx" #include "MEDCouplingCMesh.hxx" @@ -14,13 +13,18 @@ #include #include #include "vtkDataSet.h" +#include "vtkCellData.h" +#include +#include +#include using namespace MEDCoupling; void Visualization::CatalystInitialize(const std::string& script) { + if(Processor == NULL) { Processor = vtkCPProcessor::New(); @@ -46,13 +50,7 @@ void Visualization::CatalystFinalize() Processor->Delete(); Processor = NULL; } -/* uuuuuuuh ... maybe cleaning vtkDataSet -if(VTKGrid) - { - VTKGrid->Delete(); - VTKGrid = NULL; - } -*/ + return; } @@ -75,9 +73,9 @@ void Visualization::CatalystCoProcess(vtkDataSet *VTKGrid, double time, dataDescription->Delete(); } -void Visualization::ConvertToVTK(MEDCoupling::MEDCouplingFieldDouble* field, vtkDataSet *VTKGrid) +void Visualization::ConvertToVTK(MEDCoupling::MEDCouplingFieldDouble* field, vtkDataSet *&VTKGrid) +//vtkDataSet * Visualization::ConvertToVTK(MEDCoupling::MEDCouplingFieldDouble* field) { - MEDCoupling::MEDCouplingFieldDouble *f = field; MEDCoupling::MCAuto ff(MEDFileField1TS::New()); ff->setFieldNoProfileSBT(f); @@ -92,6 +90,7 @@ void Visualization::ConvertToVTK(MEDCoupling::MEDCouplingFieldDouble* field, vtk { MCAuto mmu(MEDFileUMesh::New()); mmu->setMeshAtLevel(0,dynamic_cast(m)); + mmu->forceComputationOfParts(); mm=DynamicCast(mmu); } else if(dynamic_cast(m)) @@ -108,16 +107,28 @@ void Visualization::ConvertToVTK(MEDCoupling::MEDCouplingFieldDouble* field, vtk //mm=DynamicCast(mmc); mm=0; } - throw ; + else + { + throw ; + } MCAuto ms(MEDFileMeshes::New()); ms->pushMesh(mm); + ms->getMeshesNames(); + // + int proc_id; + MPI_Comm_rank(MPI_COMM_WORLD,&proc_id); + // MEDFileFieldRepresentationTree Tree; Tree.loadInMemory(fs,ms); + + Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(0,true); + Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(1,false); + Tree.activateTheFirst();//"TS0/Fluid domain/ComSup0/TempC@@][@@P0" std::string meshName; TimeKeeper TK(0); int tmp1,tmp2; double tmp3(f->getTime(tmp1,tmp2)); - vtkDataSet *ret(Tree.buildVTKInstance(true,tmp3,meshName,TK)); + vtkDataSet *ret(Tree.buildVTKInstance(false,tmp3,meshName,TK)); VTKGrid = ret; } @@ -133,6 +144,7 @@ void Visualization::run(MEDCoupling::MEDCouplingFieldDouble* field, const std::s { int proc_id; MPI_Comm_rank(MPI_COMM_WORLD,&proc_id); + if( field == NULL) { std::cerr << "Description n° " << proc_id << ": NULL pointer" << std::endl; @@ -141,20 +153,14 @@ void Visualization::run(MEDCoupling::MEDCouplingFieldDouble* field, const std::s << field->getDescription() << std::endl; - CatalystInitialize(pathPipeline); - vtkDataSet *VTKGrid = 0; ConvertToVTK(field, VTKGrid); - // Time information for CatalystCoProcess - int iteration = 0; - int order = 0; - double time = field->getTime(iteration, order); - - CatalystCoProcess(VTKGrid, time, iteration); - + const char *fileName = pathPipeline.c_str(); + CatalystInitialize(fileName); + CatalystCoProcess(VTKGrid, 0.0, 0); CatalystFinalize(); - + return ; } diff --git a/src/Insitu/visu.hxx b/src/Insitu/visu.hxx old mode 100644 new mode 100755 index c081c2d3..30d408eb --- a/src/Insitu/visu.hxx +++ b/src/Insitu/visu.hxx @@ -16,11 +16,12 @@ class Visualization { vtkCPProcessor* Processor; - private : + //private : + public : void CatalystInitialize(const std::string& pipeline); void CatalystFinalize(); void CatalystCoProcess(vtkDataSet *VTKGrid, double time, unsigned int timeStep); - void ConvertToVTK(MEDCoupling::MEDCouplingFieldDouble* field, vtkDataSet *VTKGrid); + void ConvertToVTK(MEDCoupling::MEDCouplingFieldDouble* field, vtkDataSet *&VTKGrid); public : Visualization(); void run(MEDCoupling::MEDCouplingFieldDouble*, const std::string& pathPipeline);