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