]> SALOME platform Git repositories - modules/paravis.git/blob - src/Insitu/visu.cxx
Salome HOME
Create insitu library.
[modules/paravis.git] / src / Insitu / visu.cxx
1 #include "visu.hxx"
2 // #include "MEDCouplingFieldDouble.hxx"
3 #include "MEDFileField.hxx"
4 #include "MEDCouplingUMesh.hxx"
5 #include "MEDCouplingCMesh.hxx"
6 #include "MEDCouplingCurveLinearMesh.hxx"
7 #include "MEDFileMesh.hxx"
8 #include "MEDFileFieldRepresentationTree.hxx"
9 #include <iostream>
10 #include "mpi.h" 
11
12 #include <vtkCPDataDescription.h>
13 #include <vtkCPInputDataDescription.h>
14 #include <vtkCPPythonScriptPipeline.h>
15 #include <vtkCPProcessor.h>
16 #include "vtkDataSet.h"
17
18
19 using namespace MEDCoupling;
20
21
22 void Visualization::CatalystInitialize(const std::string& script)
23 {
24   if(Processor == NULL)
25     {
26     Processor = vtkCPProcessor::New();
27     Processor->Initialize();
28     }
29   else
30     {
31     Processor->RemoveAllPipelines();
32     }
33   // Add in the Python script
34   vtkCPPythonScriptPipeline* pipeline = vtkCPPythonScriptPipeline::New();
35   char *c = const_cast<char*>(script.c_str());
36   pipeline->Initialize(c);
37   Processor->AddPipeline(pipeline);
38   pipeline->Delete();
39   return;
40 }
41
42 void Visualization::CatalystFinalize()
43 {
44   if(Processor)
45     {
46     Processor->Delete();
47     Processor = NULL;
48     }
49 /* uuuuuuuh ... maybe cleaning vtkDataSet  
50 if(VTKGrid)
51     {
52     VTKGrid->Delete();
53     VTKGrid = NULL;
54     }
55 */
56   return;
57 }
58
59 void Visualization::CatalystCoProcess(vtkDataSet *VTKGrid, double time,
60                                       unsigned int timeStep)
61 {
62   vtkCPDataDescription* dataDescription = vtkCPDataDescription::New();
63   // specify the simulation time and time step for Catalyst
64   dataDescription->AddInput("input");
65   dataDescription->SetTimeData(time, timeStep);
66
67   if(Processor->RequestDataDescription(dataDescription) != 0)
68     {
69     // Make a map from input to our VTK grid so that
70     // Catalyst gets the proper input dataset for the pipeline.
71     dataDescription->GetInputDescriptionByName("input")->SetGrid(VTKGrid);
72     // Call Catalyst to execute the desired pipelines.
73     Processor->CoProcess(dataDescription);
74     }
75   dataDescription->Delete();
76 }
77
78 void Visualization::ConvertToVTK(MEDCoupling::MEDCouplingFieldDouble* field, vtkDataSet *VTKGrid)
79 {
80
81   MEDCoupling::MEDCouplingFieldDouble *f = field;
82   MEDCoupling::MCAuto<MEDFileField1TS> ff(MEDFileField1TS::New());
83   ff->setFieldNoProfileSBT(f);
84   MCAuto<MEDFileFieldMultiTS> fmts(MEDFileFieldMultiTS::New());
85   fmts->pushBackTimeStep(ff);
86   MCAuto<MEDFileFields> fs(MEDFileFields::New());
87   fs->pushField(fmts);
88   
89   MEDCouplingMesh *m(f->getMesh());
90   MCAuto<MEDFileMesh> mm;
91   if(dynamic_cast<MEDCouplingUMesh *>(m))
92     {
93       MCAuto<MEDFileUMesh> mmu(MEDFileUMesh::New()); 
94       mmu->setMeshAtLevel(0,dynamic_cast<MEDCouplingUMesh *>(m));
95       mm=DynamicCast<MEDFileUMesh,MEDFileMesh>(mmu);
96     }
97   else if(dynamic_cast<MEDCouplingCMesh *>(m))
98     {
99       MCAuto<MEDFileCMesh> mmc(MEDFileCMesh::New()); 
100       mmc->setMesh(dynamic_cast<MEDCouplingCMesh *>(m));
101       mm=DynamicCast<MEDFileCMesh,MEDFileMesh>(mmc);
102     }
103   else if(dynamic_cast<MEDCouplingCurveLinearMesh *>(m))
104     {
105       // MCAuto<MEDFileCLMesh> mmc(MEDFileCLMesh::New()); 
106       MCAuto<MEDFileCurveLinearMesh> mmc(MEDFileCurveLinearMesh::New()); 
107       mmc->setMesh(dynamic_cast<MEDCouplingCurveLinearMesh *>(m));
108       //mm=DynamicCast<MEDCouplingCurveLinearMesh,MEDFileMesh>(mmc);
109       mm=0;
110     }
111   throw ;
112   MCAuto<MEDFileMeshes> ms(MEDFileMeshes::New());
113   ms->pushMesh(mm);
114   MEDFileFieldRepresentationTree Tree;
115   Tree.loadInMemory(fs,ms);
116   std::string meshName;
117   TimeKeeper TK(0);
118   int tmp1,tmp2;
119   double tmp3(f->getTime(tmp1,tmp2));
120   vtkDataSet *ret(Tree.buildVTKInstance(true,tmp3,meshName,TK));
121   VTKGrid = ret;
122 }
123
124 Visualization::Visualization()
125 {
126 Processor = NULL;
127 }
128
129
130 // ajouter en parametre le fichier python qui contient le pipeline
131 // enlever le const s'il gene
132 void Visualization::run(MEDCoupling::MEDCouplingFieldDouble* field, const std::string& pathPipeline)
133 {
134   int proc_id;
135   MPI_Comm_rank(MPI_COMM_WORLD,&proc_id);
136   if( field == NULL)
137   {
138     std::cerr << "Description n° " << proc_id << ": NULL pointer" << std::endl;
139   }
140   std::cerr << "Description n° " << proc_id << ":"
141             << field->getDescription() << std::endl;
142
143
144   CatalystInitialize(pathPipeline);
145
146   vtkDataSet *VTKGrid = 0;
147   ConvertToVTK(field, VTKGrid);
148
149   // Time information for CatalystCoProcess
150   int iteration = 0;
151   int order = 0;
152   double time = field->getTime(iteration, order);
153
154   CatalystCoProcess(VTKGrid, time, iteration);
155   
156   CatalystFinalize();
157
158   return ;
159 }
160