Salome HOME
Copyright update 2022
[modules/paravis.git] / src / Plugins / ParaMEDCorba / plugin / ParaMEDMEM2VTK / VTKMEDCouplingFieldClient.cxx
1 // Copyright (C) 2010-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 "VTKMEDCouplingFieldClient.hxx"
21 #include "VTKMEDCouplingMeshClient.hxx"
22
23 #include "vtkErrorCode.h"
24 #include "vtkCellData.h"
25 #include "vtkPointData.h"
26 #include "vtkDoubleArray.h"
27 #include "vtkUnstructuredGrid.h"
28
29 #include <string>
30 #include <sstream>
31
32 std::vector<double> ParaMEDMEM2VTK::FillMEDCouplingFieldDoubleInstanceFrom(SALOME_MED::MEDCouplingFieldDoubleCorbaInterface_ptr fieldPtr, vtkDataSet *ret)
33 {
34   fieldPtr->Register();
35   //
36   SALOME_MED::MEDCouplingMeshCorbaInterface_var meshPtr=fieldPtr->getMesh();
37   ParaMEDMEM2VTK::FillMEDCouplingMeshInstanceFrom(meshPtr,ret);
38   meshPtr->UnRegister();
39   //
40   std::vector<double> ret2=FillMEDCouplingFieldDoublePartOnly(fieldPtr,ret);
41   fieldPtr->UnRegister();
42   //
43   return ret2;
44 }
45
46 std::vector<double> ParaMEDMEM2VTK::FillMEDCouplingFieldDoublePartOnly(SALOME_MED::MEDCouplingFieldDoubleCorbaInterface_ptr fieldPtr, vtkDataSet *ret)
47 {
48   std::vector<double> ret2;
49   //
50   SALOME_TYPES::ListOfLong *tinyL;
51   SALOME_TYPES::ListOfDouble *tinyD;
52   SALOME_TYPES::ListOfString *tinyS;
53   fieldPtr->getTinyInfo(tinyL,tinyD,tinyS);
54   //
55   int typeOfField=(*tinyL)[0];// 0:ON_CELLS, 1:ON_NODES, 2:ON_GAUSS_PT, 3:ON_GAUSS_NE
56   int timeDiscr=(*tinyL)[1];
57   int nbOfTuples=(*tinyL)[3];
58   int nbOfComponents=(*tinyL)[4];
59   std::vector<std::string> comps(nbOfComponents);
60   for(int i=0;i<nbOfComponents;i++)
61     {
62       std::string comp((*tinyS)[i]);
63       if(!comp.empty())
64         comps[i]=comp;
65       else
66         {
67           std::ostringstream oss; oss << "Comp" << i;
68           comps[i]=oss.str();
69         }
70     }
71   std::string name;
72   if(timeDiscr!=7)
73     name=((*tinyS)[nbOfComponents]);
74   else
75     name=((*tinyS)[2*nbOfComponents]);
76   //
77   switch(timeDiscr)
78     {
79     case 4://NO_TIME
80       ret2.resize(1);
81       ret2[0]=-1;
82       break;
83     case 5://ONE_TIME
84       ret2.resize(1);
85       ret2[0]=(*tinyD)[1];
86       break;
87     case 6://LINEAR_TIME
88     case 7://CONST_ON_TIME_INTERVAL
89       ret2.resize(1);
90       ret2[0]=((*tinyD)[1]+(*tinyD)[2])/2.;
91       break;
92     default:
93       vtkErrorWithObjectMacro(ret,"Unrecognized time discretization of Field ! Not implemented yet !");
94     }
95   //
96   delete tinyL;
97   delete tinyD;
98   delete tinyS;
99   //
100   vtkDoubleArray *var0=vtkDoubleArray::New();
101   vtkDoubleArray *var1=0;
102   double *var0Ptr=0;
103   double *var1Ptr=0;
104   var0->SetName(name.c_str());
105   var0->SetNumberOfComponents(nbOfComponents);
106   var0->SetNumberOfTuples(nbOfTuples);
107   for(int i=0;i<nbOfComponents;i++)
108     var0->SetComponentName(i,comps[i].c_str());
109   var0Ptr=var0->GetPointer(0);
110   if(timeDiscr==7)
111     {
112       var1->SetNumberOfComponents(nbOfComponents);
113       var1->SetNumberOfTuples(nbOfTuples);
114       var1Ptr=var1->GetPointer(0);
115       std::ostringstream oss; oss << name << "_end_array";
116       var1->SetName(oss.str().c_str());
117     }
118   //
119   SALOME_TYPES::ListOfLong *bigDataI;
120   SALOME_TYPES::ListOfDouble2 *bigArr;
121   fieldPtr->getSerialisationData(bigDataI,bigArr);
122   delete bigDataI;//for the moment gauss points are not dealt
123   SALOME_TYPES::ListOfDouble& oneArr=(*bigArr)[0];
124   int nbVals=oneArr.length();
125   for(int i=0;i<nbVals;i++)
126     var0Ptr[i]=oneArr[i];
127   if(timeDiscr==7)
128     {
129       SALOME_TYPES::ListOfDouble& scndArr=(*bigArr)[1];
130       for(int i=0;i<nbVals;i++)
131         var1Ptr[i]=scndArr[i];
132     }
133   delete bigArr;
134   //
135   switch(typeOfField)
136     {
137     case 0://ON_CELLS
138       ret->GetCellData()->AddArray(var0);
139       if(timeDiscr==7)
140         ret->GetCellData()->AddArray(var1);
141       break;
142     case 1://ON_NODES
143       ret->GetPointData()->AddArray(var0);
144       if(timeDiscr==7)
145         ret->GetCellData()->AddArray(var1);
146       break;
147     default:
148       vtkErrorWithObjectMacro(ret,"Not implemented yet for gauss fields and gauss on elements fields !");
149     }
150   var0->Delete();
151   if(timeDiscr==7)
152     var1->Delete();
153   //
154   return ret2;
155 }
156
157 vtkDataSet *ParaMEDMEM2VTK::BuildFullyFilledFromMEDCouplingFieldDoubleInstance(SALOME_MED::MEDCouplingFieldDoubleCorbaInterface_ptr fieldPtr,
158                                                                                std::vector<double>& times)
159 {
160   fieldPtr->Register();
161   SALOME_MED::MEDCouplingMeshCorbaInterface_var meshPtr=fieldPtr->getMesh();
162   bool dummy;//VTK bug
163   vtkDataSet *ret=ParaMEDMEM2VTK::BuildFromMEDCouplingMeshInstance(meshPtr,dummy);//VTK bug
164   meshPtr->UnRegister();
165   //
166   std::vector<double> ret2=FillMEDCouplingFieldDoublePartOnly(fieldPtr,ret);
167   times=ret2;
168   //
169   fieldPtr->UnRegister();
170   return ret;
171 }
172
173 vtkDoubleArray *ParaMEDMEM2VTK::BuildFromMEDCouplingFieldDoubleArr(SALOME_MED::DataArrayDoubleCorbaInterface_ptr dadPtr)
174 {
175   vtkDoubleArray *ret=vtkDoubleArray::New();
176   //
177   SALOME_TYPES::ListOfLong *tinyL=0;
178   SALOME_TYPES::ListOfDouble *bigD=0;
179   SALOME_TYPES::ListOfString *tinyS=0;
180   //
181   dadPtr->getTinyInfo(tinyL,tinyS);
182   int nbOfTuples=(*tinyL)[0];
183   int nbOfCompo=(*tinyL)[1];
184   std::string name((*tinyS)[0]);
185   std::vector<std::string> comps(nbOfCompo);
186   for(int i=0;i<nbOfCompo;i++)
187     comps[i]=(*tinyS)[i+1];
188   delete tinyL;
189   delete tinyS;
190   //
191   ret->SetName(name.c_str());
192   ret->SetNumberOfComponents(nbOfCompo);
193   ret->SetNumberOfTuples(nbOfTuples);
194   for(int i=0;i<nbOfCompo;i++)
195     {
196       if(!comps[i].empty())
197         ret->SetComponentName(i,comps[i].c_str());
198       else
199         {
200           std::ostringstream oss; oss << "Comp" << i;
201           ret->SetComponentName(i,oss.str().c_str());
202         }
203     }
204   int nbElems=nbOfCompo*nbOfTuples;
205   double *pt=ret->GetPointer(0);
206   dadPtr->getSerialisationData(bigD);
207   for(int i=0;i<nbElems;i++)
208     pt[i]=(*bigD)[i];
209   delete bigD;
210   //
211   return ret;
212 }