Salome HOME
Copyright update 2022
[modules/paravis.git] / src / Plugins / ParaMEDCorba / plugin / ParaMEDMEM2VTK / VTKMEDCouplingMultiFieldsClient.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 // To access to vtkUnstructuredGrid::Faces and FaceLocations
21 #ifdef WIN32
22 #define _ALLOW_KEYWORD_MACROS
23 #endif
24 #define protected public
25
26 #include "VTKMEDCouplingMultiFieldsClient.hxx"
27 #include "VTKMEDCouplingMeshClient.hxx"
28 #include "VTKMEDCouplingFieldClient.hxx"
29
30 #include "vtkUnstructuredGrid.h"
31 #include "vtkRectilinearGrid.h"
32 #include "vtkDoubleArray.h"
33 #include "vtkErrorCode.h"
34 #include "vtkCellData.h"
35 #include "vtkIdTypeArray.h"
36 #include "vtkPointData.h"
37
38 #include <sstream>
39 #include <iterator>
40 #include <algorithm>
41 #include <functional>
42
43 const double ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::EPS_TIME=1e-7;
44
45 ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::MEDCouplingMultiFieldsFetcher(int bufferingPolicy,
46     SALOME_MED::MEDCouplingMultiFieldsCorbaInterface_ptr mfieldsPtr):_effective_pol(bufferingPolicy),_mfields_ptr_released(false)
47 {
48   _mfields_ptr=SALOME_MED::MEDCouplingMultiFieldsCorbaInterface::_duplicate(mfieldsPtr);
49   _mfields_ptr->Register();
50 }
51
52 ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::~MEDCouplingMultiFieldsFetcher()
53 {
54   for(std::vector<vtkDataSet *>::iterator it=_meshes.begin();it!=_meshes.end();it++)
55   {
56     if(*it)
57       (*it)->Delete();
58   }
59   for(std::vector<vtkDoubleArray *>::iterator it2=_arrays.begin();it2!=_arrays.end();it2++)
60   {
61     if(*it2)
62       (*it2)->Delete();
63   }
64   if(!_mfields_ptr_released)
65     _mfields_ptr->UnRegister();
66 }
67
68 std::vector<double> ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::getTimeStepsForPV()
69 {
70   retrievesMainTinyInfo();
71   int nbOfFields=_mesh_id_per_field.size();
72   //
73   _time_label_per_field.resize(nbOfFields);
74   SALOME_MED::MEDCouplingFieldOverTimeCorbaInterface_var fotPtr=SALOME_MED::MEDCouplingFieldOverTimeCorbaInterface::_narrow(_mfields_ptr);
75   if(CORBA::is_nil(fotPtr))
76   {
77     for(CORBA::Long i=0;i<nbOfFields;i++)
78       _time_label_per_field[i]=(double)i;
79   }
80   else
81   {
82     double tmp=0.;
83     for(CORBA::Long i=0;i<nbOfFields;i++)
84     {
85       if(!_time_def_per_field[i].empty())
86         _time_label_per_field[i]=_time_def_per_field[i].front();
87       else
88         _time_label_per_field[i]=tmp++;
89     }
90   }
91   return _time_label_per_field;
92 }
93
94 void ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::fetchRegardingPolicy()
95 {
96   if(_effective_pol>=10)
97   {
98     fetchAll();
99     return ;
100   }
101   if(_effective_pol>=1 && _effective_pol<=9)
102   {
103     fetchMeshes();
104     return ;
105   }
106 }
107
108 vtkDataSet *ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::buildDataSetOnTime(double time)
109 {
110   int fieldId=getPosGivenTimeLabel(time);
111   if(fieldId<0)
112     return 0;
113   fetchDataIfNeeded(fieldId);
114   int meshId=_mesh_id_per_field[fieldId];
115   vtkDataSet *ret0=_meshes[meshId];
116   std::string clsName=ret0->GetClassName();
117   if(clsName=="vtkUnstructuredGrid")
118   {
119     vtkUnstructuredGrid *ret1=vtkUnstructuredGrid::New();
120     ret1->DeepCopy(ret0);
121     if(_is_meshes_polyhedron[meshId])//bug VTK polyhedron
122     {//bug VTK polyhedron part
123       ret1->Faces->UnRegister(ret1);
124       ret1->Faces=vtkIdTypeArray::New();
125       ret1->Faces->DeepCopy(((vtkUnstructuredGrid *)ret0)->GetFaces());
126       ret1->Faces->Register(ret1);
127       ret1->Faces->Delete();
128       ret1->FaceLocations->UnRegister(ret1);
129       ret1->FaceLocations=vtkIdTypeArray::New();
130       ret1->FaceLocations->DeepCopy(((vtkUnstructuredGrid *)ret0)->GetFaceLocations());
131       ret1->FaceLocations->Register(ret1);
132       ret1->FaceLocations->Delete();
133     }//end bug VTK polyhedron part
134     appendFieldValueOnAlreadyFetchedData(ret1,fieldId);
135     applyBufferingPolicy();
136     return ret1;
137   }
138   if(clsName=="vtkRectilinearGrid")
139   {
140     vtkRectilinearGrid *ret1=vtkRectilinearGrid::New();
141     ret1->DeepCopy(ret0);
142     appendFieldValueOnAlreadyFetchedData(ret1,fieldId);
143     applyBufferingPolicy();
144     return ret1;
145   }
146   return 0;
147 }
148
149 void ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::retrievesMainTinyInfo()
150 {
151   SALOME_TYPES::ListOfLong *tinyL=0;
152   SALOME_TYPES::ListOfDouble *tinyD=0;
153   SALOME_TYPES::ListOfString *tinyS=0;
154   //
155   CORBA::Long nbOfArrays;
156   CORBA::Long nbOfFields;
157   CORBA::Long nbOfMeshes=_mfields_ptr->getMainTinyInfo(tinyL,tinyD,nbOfArrays,nbOfFields);
158   int sz=(*tinyL)[0];//nbOfFields
159   //int sz2=(*tinyL)[1];//sigma(nbOfArraysPerField)
160   _time_discr_per_field.resize(sz);//4 : NO_TIME  5:ONE_TIME 6:LINEAR_TIME 7:CONST_ON_TIME_INTERVAL
161   _mesh_id_per_field.resize(sz);
162   _array_ids_per_field.resize(sz);
163   _time_def_per_field.resize(sz);
164   int offsetTime=0;
165   int offsetArrays=0;
166   for(int i=0;i<sz;i++)
167   {
168     _mesh_id_per_field[i]=(*tinyL)[3+i];
169     int nbOfArrayForCurField=(*tinyL)[sz+3+i];
170     _array_ids_per_field[i].resize(nbOfArrayForCurField);
171     for(int k=0;k<nbOfArrayForCurField;k++)
172       _array_ids_per_field[i][k]=(*tinyL)[5*sz+3+offsetArrays+k];
173     _time_discr_per_field[i]=(*tinyL)[2*sz+3+i];
174     int nbOfTimeSpot=(*tinyL)[3*sz+3+i]-1;//-1 because time precision is not useful here.
175     _time_def_per_field[i].resize(nbOfTimeSpot);
176     for(int j=0;j<nbOfTimeSpot;j++)
177       _time_def_per_field[i][j]=(*tinyD)[offsetTime+1+j];
178     offsetTime+=nbOfTimeSpot+1;
179     offsetArrays+=nbOfArrayForCurField;
180   }
181   delete tinyL;
182   delete tinyD;
183   //
184   _meshes.resize(nbOfMeshes+1);
185   _is_meshes_polyhedron.resize(nbOfMeshes+1);
186   _arrays.resize(nbOfArrays+1);
187   //
188   _info_per_field.resize(nbOfFields);
189   for(int i=0;i<nbOfFields;i++)
190   {
191     _mfields_ptr->getTinyInfo(i,tinyL,tinyD,tinyS);
192     _info_per_field[i]._type=(*tinyL)[0];
193     _info_per_field[i]._name=(*tinyS)[0];
194     delete tinyL;
195     delete tinyD;
196     delete tinyS;
197   }
198 }
199
200 void ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::fetchAll()
201 {
202   fetchMeshes();
203   int nbOfArrays=_arrays.size();
204   for(int i=0;i<nbOfArrays;i++)
205   {
206     SALOME_MED::DataArrayDoubleCorbaInterface_var daPtr=_mfields_ptr->getArray(i);
207     if(_arrays[i])
208       _arrays[i]->Delete();
209     _arrays[i]=ParaMEDMEM2VTK::BuildFromMEDCouplingFieldDoubleArr(daPtr);
210     daPtr->UnRegister();
211   }
212   unregisterRemoteServantIfAllFetched();
213 }
214
215 /*!
216  * Fetches meshes without regarding if already fetched
217  */
218 void ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::fetchMeshes()
219 {
220   int nbOfMeshes=_meshes.size();
221   for(int i=0;i<nbOfMeshes;i++)
222   {
223     SALOME_MED::MEDCouplingMeshCorbaInterface_var mPtr=_mfields_ptr->getMeshWithId(i);
224     if(_meshes[i])
225       _meshes[i]->Delete();
226     bool polyh=false;//bug VTK
227     _meshes[i]=ParaMEDMEM2VTK::BuildFromMEDCouplingMeshInstance(mPtr,polyh);//bug VTK
228     _is_meshes_polyhedron[i]=polyh;//bug VTK
229     mPtr->UnRegister();
230   }
231   unregisterRemoteServantIfAllFetched();
232 }
233
234 /*!
235  * For a field with id 'fieldId' this method CORBA fetch, if needed, basic data.
236  * 'fieldId' should be correct no check of that is done !
237  */
238 void ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::fetchDataIfNeeded(int fieldId)
239 {
240   std::vector<int> arrayIds=_array_ids_per_field[fieldId];
241   int meshId=_mesh_id_per_field[fieldId];
242   if(!_meshes[meshId])
243   {
244     SALOME_MED::MEDCouplingMeshCorbaInterface_var mPtr=_mfields_ptr->getMeshWithId(meshId);
245     bool polyh=false;//bug VTK
246     _meshes[meshId]=ParaMEDMEM2VTK::BuildFromMEDCouplingMeshInstance(mPtr,polyh);//bug VTK
247     _is_meshes_polyhedron[meshId]=polyh;//bug VTK
248     mPtr->UnRegister();
249   }
250   for(std::vector<int>::const_iterator it=arrayIds.begin();it!=arrayIds.end();it++)
251   {
252     if(!_arrays[*it])
253     {
254       SALOME_MED::DataArrayDoubleCorbaInterface_var daPtr=_mfields_ptr->getArray(*it);
255       _arrays[*it]=ParaMEDMEM2VTK::BuildFromMEDCouplingFieldDoubleArr(daPtr);
256       daPtr->UnRegister();
257     }
258   }
259   unregisterRemoteServantIfAllFetched();
260 }
261
262 void ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::unregisterRemoteServantIfAllFetched()
263 {
264   for(std::vector<vtkDataSet *>::iterator it=_meshes.begin();it!=_meshes.end();it++)
265   {
266     if((*it)==0)
267       return ;
268   }
269   for(std::vector<vtkDoubleArray *>::iterator it2=_arrays.begin();it2!=_arrays.end();it2++)
270   {
271     if((*it2)==0)
272       return ;
273   }
274   if(!_mfields_ptr_released)
275   {
276     _mfields_ptr_released=true;
277     _mfields_ptr->UnRegister();
278   }
279 }
280
281 void ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::applyBufferingPolicy()
282 {
283   if(_effective_pol==0)
284   {//
285     for(std::vector<vtkDataSet *>::iterator it=_meshes.begin();it!=_meshes.end();it++)
286     {
287       if(*it)
288       {
289         (*it)->Delete();
290         *it=0;
291       }
292     }
293     for(std::vector<vtkDoubleArray *>::iterator it2=_arrays.begin();it2!=_arrays.end();it2++)
294     {
295       if(*it2)
296       {
297         (*it2)->Delete();
298         *it2=0;
299       }
300     }
301   }
302   //else nothing to do let the plugin bufferize
303 }
304
305 void ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::appendFieldValueOnAlreadyFetchedData(vtkDataSet *ds, int fieldId)
306 {
307   const TinyInfoOnField& info=_info_per_field[fieldId];
308   vtkDoubleArray *arr=_arrays[_array_ids_per_field[fieldId].front()];
309   arr->SetName(info._name.c_str());
310   if(info._type==0)//ON_CELLS
311   {
312     ds->GetCellData()->AddArray(arr);
313     return ;
314   }
315   if(info._type==1)//ON_NODES
316   {
317     ds->GetPointData()->AddArray(arr);
318     return ;
319   }
320 }
321
322 int ParaMEDMEM2VTK::MEDCouplingMultiFieldsFetcher::getPosGivenTimeLabel(double t)
323 {
324   int nbOfFields=_time_label_per_field.size();
325   for(int i=0;i<nbOfFields;i++)
326     if(fabs(_time_label_per_field[i]-t)<EPS_TIME)
327       return i;
328   //2nd chance
329   std::vector<double>::iterator it=std::find_if(_time_label_per_field.begin(),_time_label_per_field.end(),
330       std::bind(std::greater<double>(),std::placeholders::_1,t));
331   if(it!=_time_label_per_field.end() && it!=_time_label_per_field.end())
332     return std::distance(_time_label_per_field.begin(),it);
333   //
334   std::ostringstream oss;
335   oss << "Unexisting time : " << t << " Not in ";
336   std::copy(_time_label_per_field.begin(),_time_label_per_field.end(),std::ostream_iterator<double>(oss," "));
337   oss << " !";
338   vtkOutputWindowDisplayErrorText(oss.str().c_str());
339   return -1;
340 }