1 // Copyright (C) 2010-2016 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
20 #include "VTKMEDCouplingFieldClient.hxx"
21 #include "VTKMEDCouplingMeshClient.hxx"
23 #include "vtkErrorCode.h"
24 #include "vtkCellData.h"
25 #include "vtkPointData.h"
26 #include "vtkDoubleArray.h"
27 #include "vtkUnstructuredGrid.h"
32 std::vector<double> ParaMEDMEM2VTK::FillMEDCouplingFieldDoubleInstanceFrom(SALOME_MED::MEDCouplingFieldDoubleCorbaInterface_ptr fieldPtr, vtkDataSet *ret)
36 SALOME_MED::MEDCouplingMeshCorbaInterface_var meshPtr=fieldPtr->getMesh();
37 ParaMEDMEM2VTK::FillMEDCouplingMeshInstanceFrom(meshPtr,ret);
38 meshPtr->UnRegister();
40 std::vector<double> ret2=FillMEDCouplingFieldDoublePartOnly(fieldPtr,ret);
41 fieldPtr->UnRegister();
46 std::vector<double> ParaMEDMEM2VTK::FillMEDCouplingFieldDoublePartOnly(SALOME_MED::MEDCouplingFieldDoubleCorbaInterface_ptr fieldPtr, vtkDataSet *ret)
48 std::vector<double> ret2;
50 SALOME_TYPES::ListOfLong *tinyL;
51 SALOME_TYPES::ListOfDouble *tinyD;
52 SALOME_TYPES::ListOfString *tinyS;
53 fieldPtr->getTinyInfo(tinyL,tinyD,tinyS);
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++)
62 std::string comp((*tinyS)[i]);
67 std::ostringstream oss; oss << "Comp" << i;
73 name=((*tinyS)[nbOfComponents]);
75 name=((*tinyS)[2*nbOfComponents]);
88 case 7://CONST_ON_TIME_INTERVAL
90 ret2[0]=((*tinyD)[1]+(*tinyD)[2])/2.;
93 vtkErrorWithObjectMacro(ret,"Unrecognized time discretization of Field ! Not implemented yet !");
100 vtkDoubleArray *var0=vtkDoubleArray::New();
101 vtkDoubleArray *var1=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);
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());
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];
129 SALOME_TYPES::ListOfDouble& scndArr=(*bigArr)[1];
130 for(int i=0;i<nbVals;i++)
131 var1Ptr[i]=scndArr[i];
138 ret->GetCellData()->AddArray(var0);
140 ret->GetCellData()->AddArray(var1);
143 ret->GetPointData()->AddArray(var0);
145 ret->GetCellData()->AddArray(var1);
148 vtkErrorWithObjectMacro(ret,"Not implemented yet for gauss fields and gauss on elements fields !");
157 vtkDataSet *ParaMEDMEM2VTK::BuildFullyFilledFromMEDCouplingFieldDoubleInstance(SALOME_MED::MEDCouplingFieldDoubleCorbaInterface_ptr fieldPtr,
158 std::vector<double>& times)
160 fieldPtr->Register();
161 SALOME_MED::MEDCouplingMeshCorbaInterface_var meshPtr=fieldPtr->getMesh();
163 vtkDataSet *ret=ParaMEDMEM2VTK::BuildFromMEDCouplingMeshInstance(meshPtr,dummy);//VTK bug
164 meshPtr->UnRegister();
166 std::vector<double> ret2=FillMEDCouplingFieldDoublePartOnly(fieldPtr,ret);
169 fieldPtr->UnRegister();
173 vtkDoubleArray *ParaMEDMEM2VTK::BuildFromMEDCouplingFieldDoubleArr(SALOME_MED::DataArrayDoubleCorbaInterface_ptr dadPtr)
175 vtkDoubleArray *ret=vtkDoubleArray::New();
177 SALOME_TYPES::ListOfLong *tinyL=0;
178 SALOME_TYPES::ListOfDouble *bigD=0;
179 SALOME_TYPES::ListOfString *tinyS=0;
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];
191 ret->SetName(name.c_str());
192 ret->SetNumberOfComponents(nbOfCompo);
193 ret->SetNumberOfTuples(nbOfTuples);
194 for(int i=0;i<nbOfCompo;i++)
196 if(!comps[i].empty())
197 ret->SetComponentName(i,comps[i].c_str());
200 std::ostringstream oss; oss << "Comp" << i;
201 ret->SetComponentName(i,oss.str().c_str());
204 int nbElems=nbOfCompo*nbOfTuples;
205 double *pt=ret->GetPointer(0);
206 dadPtr->getSerialisationData(bigD);
207 for(int i=0;i<nbElems;i++)