Salome HOME
Copyright update 2022
[modules/paravis.git] / src / Insitu / VisualizationLibrary / visu.cxx
1 // Copyright (C) 2016-2022  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, or (at your option) any later version.
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
20 #include "visu.hxx"
21 #include "MEDFileField.hxx"
22 #include "MEDCouplingUMesh.hxx"
23 #include "MEDCouplingCMesh.hxx"
24 #include "MEDCouplingCurveLinearMesh.hxx"
25 #include "MEDFileMesh.hxx"
26 #include "MEDFileFieldRepresentationTree.hxx"
27 #include <iostream>
28 #include "mpi.h" 
29
30 #include <vtkCPDataDescription.h>
31 #include <vtkCPInputDataDescription.h>
32 #include <vtkCPPythonScriptPipeline.h>
33 #include <vtkCPProcessor.h>
34 #include "vtkDataSet.h"
35 #include "vtkCellData.h"
36
37 #include <vtkSmartPointer.h>
38 #include <vtkXMLUnstructuredGridWriter.h>
39 #include <vtkUnstructuredGrid.h>
40
41 using namespace MEDCoupling;
42
43
44 void Visualization::CatalystInitialize(const std::string& script)
45 {
46
47   if(Processor == NULL)
48     {
49     Processor = vtkCPProcessor::New();
50     Processor->Initialize();
51     }
52   else
53     {
54     Processor->RemoveAllPipelines();
55     }
56   // Add in the Python script
57   vtkCPPythonScriptPipeline* pipeline = vtkCPPythonScriptPipeline::New();
58   char *c = const_cast<char*>(script.c_str());
59   pipeline->Initialize(c);
60   Processor->AddPipeline(pipeline);
61   pipeline->Delete();
62   return;
63 }
64
65 void Visualization::CatalystFinalize()
66 {
67   if(Processor)
68     {
69     Processor->Delete();
70     Processor = NULL;
71     }
72
73   return;
74 }
75
76 void Visualization::CatalystCoProcess(vtkDataSet *VTKGrid, double time,
77                                       unsigned int timeStep)
78 {
79   vtkCPDataDescription* dataDescription = vtkCPDataDescription::New();
80   // specify the simulation time and time step for Catalyst
81   dataDescription->AddInput("input");
82   dataDescription->SetTimeData(time, timeStep);
83
84   if(Processor->RequestDataDescription(dataDescription) != 0)
85     {
86     // Make a map from input to our VTK grid so that
87     // Catalyst gets the proper input dataset for the pipeline.
88     dataDescription->GetInputDescriptionByName("input")->SetGrid(VTKGrid);
89     // Call Catalyst to execute the desired pipelines.
90     Processor->CoProcess(dataDescription);
91     }
92   dataDescription->Delete();
93 }
94
95 void Visualization::ConvertToVTK(MEDCoupling::MEDCouplingFieldDouble* field, vtkDataSet *&VTKGrid)
96 //vtkDataSet * Visualization::ConvertToVTK(MEDCoupling::MEDCouplingFieldDouble* field)
97 {
98   MEDCoupling::MEDCouplingFieldDouble *f = field;
99   MEDCoupling::MCAuto<MEDFileField1TS> ff(MEDFileField1TS::New());
100   ff->setFieldNoProfileSBT(f);
101   MCAuto<MEDFileFieldMultiTS> fmts(MEDFileFieldMultiTS::New());
102   fmts->pushBackTimeStep(ff);
103   MCAuto<MEDFileFields> fs(MEDFileFields::New());
104   fs->pushField(fmts);
105   
106   MEDCouplingMesh *m(f->getMesh());
107   MCAuto<MEDFileMesh> mm;
108   if(dynamic_cast<MEDCouplingUMesh *>(m))
109     {
110       MCAuto<MEDFileUMesh> mmu(MEDFileUMesh::New()); 
111       mmu->setMeshAtLevel(0,dynamic_cast<MEDCouplingUMesh *>(m));
112       mmu->forceComputationOfParts();
113       mm=DynamicCast<MEDFileUMesh,MEDFileMesh>(mmu);
114     }
115   else if(dynamic_cast<MEDCouplingCMesh *>(m))
116     {
117       MCAuto<MEDFileCMesh> mmc(MEDFileCMesh::New()); 
118       mmc->setMesh(dynamic_cast<MEDCouplingCMesh *>(m));
119       mm=DynamicCast<MEDFileCMesh,MEDFileMesh>(mmc);
120     }
121   else if(dynamic_cast<MEDCouplingCurveLinearMesh *>(m))
122     {
123       // MCAuto<MEDFileCLMesh> mmc(MEDFileCLMesh::New()); 
124       MCAuto<MEDFileCurveLinearMesh> mmc(MEDFileCurveLinearMesh::New()); 
125       mmc->setMesh(dynamic_cast<MEDCouplingCurveLinearMesh *>(m));
126       //mm=DynamicCast<MEDCouplingCurveLinearMesh,MEDFileMesh>(mmc);
127       mm=0;
128     }
129   else
130     {
131     throw ;
132     }
133   MCAuto<MEDFileMeshes> ms(MEDFileMeshes::New());
134   ms->pushMesh(mm);
135   ms->getMeshesNames();
136   //
137   int proc_id;
138   MPI_Comm_rank(MPI_COMM_WORLD,&proc_id);
139   //
140   MEDFileFieldRepresentationTree Tree;
141   Tree.loadInMemory(fs,ms);
142   
143   Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(0,true);
144   Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(1,false);
145   Tree.activateTheFirst();//"TS0/Fluid domain/ComSup0/TempC@@][@@P0"
146   std::string meshName;
147   TimeKeeper TK(0);
148   int tmp1,tmp2;
149   double tmp3(f->getTime(tmp1,tmp2));
150   vtkDataSet *ret(Tree.buildVTKInstance(false,tmp3,meshName,TK));
151   VTKGrid = ret;
152 }
153
154 Visualization::Visualization()
155 {
156 Processor = NULL;
157 }
158
159
160 // ajouter en parametre le fichier python qui contient le pipeline
161 // enlever le const s'il gene
162 void Visualization::run(MEDCoupling::MEDCouplingFieldDouble* field, const std::string& pathPipeline)
163 {
164   int proc_id;
165   MPI_Comm_rank(MPI_COMM_WORLD,&proc_id);
166
167   vtkDataSet *VTKGrid = 0;
168   ConvertToVTK(field, VTKGrid);
169
170   const char *fileName = pathPipeline.c_str();
171   CatalystInitialize(fileName);
172   CatalystCoProcess(VTKGrid, 0.0, 0);
173   CatalystFinalize();
174   
175   return ;
176 }
177