Salome HOME
404180ed0b16dd84a134308ccef34d1f24a32184
[modules/paravis.git] / src / Plugins / SimpleMode / IO / vtkSimpleMode.cxx
1 // Copyright (C) 2017  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 // Author : Anthony Geay (EDF R&D)
20
21 #include "vtkSimpleMode.h"
22
23 #include "vtkAdjacentVertexIterator.h"
24 #include "vtkDataSetSurfaceFilter.h"
25 #include "vtkIntArray.h"
26 #include "vtkCellData.h"
27 #include "vtkPointData.h"
28
29 #include "vtkStreamingDemandDrivenPipeline.h"
30 #include "vtkUnstructuredGrid.h"
31 #include "vtkPolyData.h"
32 #include  "vtkMultiBlockDataSet.h"
33
34 #include "vtkInformationStringKey.h"
35 #include "vtkAlgorithmOutput.h"
36 #include "vtkObjectFactory.h"
37 #include "vtkMutableDirectedGraph.h"
38 #include "vtkMultiBlockDataSet.h"
39 #include "vtkDataSet.h"
40 #include "vtkInformationVector.h"
41 #include "vtkInformation.h"
42 #include "vtkDataArraySelection.h"
43 #include "vtkTimeStamp.h"
44 #include "vtkInEdgeIterator.h"
45 #include "vtkInformationDataObjectKey.h"
46 #include "vtkExecutive.h"
47 #include "vtkVariantArray.h"
48 #include "vtkStringArray.h"
49 #include "vtkDoubleArray.h"
50 #include "vtkCharArray.h"
51 #include "vtkUnsignedCharArray.h"
52 #include "vtkDataSetAttributes.h"
53 #include "vtkDemandDrivenPipeline.h"
54 #include "vtkDataObjectTreeIterator.h"
55 #include "vtkWarpScalar.h"
56 #include "vtkWarpVector.h"
57 #include "vtkMultiBlockDataGroupFilter.h"
58 #include "vtkCompositeDataToUnstructuredGridFilter.h"
59
60 #include <map>
61 #include <deque>
62 #include <sstream>
63
64 vtkStandardNewMacro(vtkSimpleMode);
65
66 static const char ZE_DISPLACEMENT_NAME1[]="@@ForReal?@@";
67
68 static const char ZE_DISPLACEMENT_NAME2[]="@@ForImag?@@";
69
70 static const char ZE_DISPLACEMENT_NAME3[]="MagnitudeOfCpxDisp";
71
72 static const double EPS=1e-12;
73
74 ///////////////////
75
76 class MZCException : public std::exception
77 {
78 public:
79   MZCException(const std::string& s):_reason(s) { }
80   virtual const char *what() const throw() { return _reason.c_str(); }
81   virtual ~MZCException() throw() { }
82 private:
83   std::string _reason;
84 };
85
86 vtkSmartPointer<vtkDoubleArray> ForceTo3Compo(vtkDoubleArray *arr)
87 {
88   if(!arr)
89     return vtkSmartPointer<vtkDoubleArray>();
90   int nbCompo(arr->GetNumberOfComponents()),nbTuples(arr->GetNumberOfTuples());
91   if(nbCompo==3)
92     {
93       vtkSmartPointer<vtkDoubleArray> ret(arr);
94       arr->Register(0);
95       return ret;
96     }
97   if(nbCompo==6)
98     {
99       vtkSmartPointer<vtkDoubleArray> ret(vtkSmartPointer<vtkDoubleArray>::New());
100       ret->SetNumberOfComponents(3);
101       ret->SetNumberOfTuples(nbTuples);
102       const double *srcPt(arr->Begin());
103       double *destPt(ret->Begin());
104       for(int i=0;i<nbTuples;i++,destPt+=3,srcPt+=6)
105         std::copy(srcPt,srcPt+3,destPt);
106       return ret;
107     }
108   throw MZCException("ForceTo3Compo : internal error ! 6 or 3 compo arrays expected !");
109 }
110
111 std::vector< std::string > GetPossibleArrayNames(vtkDataSet *dataset)
112 {
113   if(!dataset)
114     throw MZCException("The input dataset is null !");
115   std::vector< std::string > ret;
116   vtkPointData *att(dataset->GetPointData());
117   for(int i=0;i<att->GetNumberOfArrays();i++)
118     {
119       vtkDataArray *locArr(att->GetArray(i));
120       int nbComp(locArr->GetNumberOfComponents());
121       if(nbComp!=3 && nbComp!=6)
122         continue;
123       std::string s(locArr->GetName());
124       ret.push_back(s);
125     }
126   return ret;
127 }
128
129 std::string FindTheBest(const std::vector<std::string>& arrNames, const std::string& key0, const std::string& key1)
130 {
131   std::string ret;
132   char points(0);
133   if(arrNames.empty())
134     return ret;
135   for(std::vector<std::string>::const_iterator it=arrNames.begin();it!=arrNames.end();it++)
136     {
137       char curNbPts(1);
138       if((*it).find(key0,0)!=std::string::npos)
139         curNbPts++;
140       if((*it).find(key1,0)!=std::string::npos)
141         curNbPts++;
142       if(curNbPts>points)
143         {
144           points=curNbPts;
145           ret=*it;
146         }
147     }
148   return ret;
149 }
150
151 std::string FindBestRealAmong(const std::vector<std::string>& arrNames)
152 {
153   static const char KEY1[]="DEPL";
154   static const char KEY2[]="REEL";
155   return FindTheBest(arrNames,KEY1,KEY2);
156 }
157
158 std::string FindBestImagAmong(const std::vector<std::string>& arrNames)
159 {
160   static const char KEY1[]="DEPL";
161   static const char KEY2[]="IMAG";
162   return FindTheBest(arrNames,KEY1,KEY2);
163 }
164
165 vtkUnstructuredGrid *ExtractInfo1(vtkInformationVector *inputVector)
166 {
167   vtkInformation *inputInfo(inputVector->GetInformationObject(0));
168   vtkDataSet *input(0);
169   vtkDataSet *input0(vtkDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
170   vtkMultiBlockDataSet *input1(vtkMultiBlockDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
171   if(input0)
172     input=input0;
173   else
174     {
175       if(!input1)
176         throw MZCException("Input dataSet must be a DataSet or single elt multi block dataset expected !");
177       if(input1->GetNumberOfBlocks()!=1)
178         throw MZCException("Input dataSet is a multiblock dataset with not exactly one block ! Use MergeBlocks or ExtractBlocks filter before calling this filter !");
179       vtkDataObject *input2(input1->GetBlock(0));
180       if(!input2)
181         throw MZCException("Input dataSet is a multiblock dataset with exactly one block but this single element is NULL !");
182       vtkDataSet *input2c(vtkDataSet::SafeDownCast(input2));
183       if(!input2c)
184         throw MZCException("Input dataSet is a multiblock dataset with exactly one block but this single element is not a dataset ! Use MergeBlocks or ExtractBlocks filter before calling this filter !");
185       input=input2c;
186     }
187   if(!input)
188     throw MZCException("Input data set is NULL !");
189   vtkUnstructuredGrid *usgIn(vtkUnstructuredGrid::SafeDownCast(input));
190   if(!usgIn)
191     throw MZCException("Input data set is not an unstructured mesh ! This filter works only on unstructured meshes !");
192   return usgIn;
193 }
194
195 void ExtractInfo3(vtkDataSet *ds, const std::string& arrName, vtkDataArray *& arr, int& idx)
196 {
197   vtkPointData *att(ds->GetPointData());
198   if(!att)
199     throw MZCException("Input dataset has no point data attribute ! Impossible to move mesh !");
200   vtkDataArray *zeArr(0);
201   int i(0);
202   for(;i<att->GetNumberOfArrays();i++)
203     {
204       vtkDataArray *locArr(att->GetArray(i));
205       std::string s(locArr->GetName());
206       if(s==arrName)
207         {
208           zeArr=locArr;
209           break;
210         }
211     }
212   if(!zeArr)
213     {
214       std::ostringstream oss;
215       oss << "Impossible to locate the array called \"" << arrName << "\" used to move mesh !";
216       throw MZCException(oss.str());
217     }
218   arr=zeArr;
219   idx=i;
220 }
221
222   
223 void ExtractInfo2(vtkDataSet *ds, const std::string& arrName, vtkDoubleArray *& arr)
224 {
225   vtkDataArray *zeArr(0);
226   int dummy;
227   ExtractInfo3(ds,arrName,zeArr,dummy);
228   arr=vtkDoubleArray::SafeDownCast(zeArr);
229   if(!arr)
230     {
231       std::ostringstream oss;
232       oss << "Array called \"" << arrName << "\" has been located but this is NOT a float64 array !";
233       throw MZCException(oss.str());
234     }
235   if(arr->GetNumberOfComponents()!=3 && arr->GetNumberOfComponents()!=6)
236     {
237       std::ostringstream oss;
238       oss << "Float64 array called \"" << arrName << "\" has been located but this array has not exactly 3 or 6 components as it should !";
239       throw MZCException(oss.str());
240     }
241   if(arr->GetNumberOfTuples()!=ds->GetNumberOfPoints())
242     {
243       std::ostringstream oss;
244       oss << "Float64-1 components array called \"" << arrName << "\" has been located but the number of tuples is invalid ! Should be " << ds->GetNumberOfPoints() << " instead of " << arr->GetNumberOfTuples() << " !";
245       throw MZCException(oss.str());
246     }
247 }
248
249 ////////////////////
250
251 class vtkSimpleMode::vtkSimpleModeInternal
252 {
253 public:
254   vtkSimpleModeInternal():_surface(0) { }
255   vtkPolyData *performConnection(vtkDataSet *ds);
256   void setFieldForReal(const std::string& st) { _real=st; }
257   std::string getFieldForReal() const { return _real; }
258   ~vtkSimpleModeInternal();
259 private:
260   vtkDataSetSurfaceFilter *_surface;
261 private:
262   std::string _real;
263 };
264
265 vtkPolyData *vtkSimpleMode::vtkSimpleModeInternal::performConnection(vtkDataSet *ds)
266 {
267   if(!_surface)
268     {
269       _surface=vtkDataSetSurfaceFilter::New();
270       _surface->SetInputData(ds);
271     }
272   _surface->Update();
273   return _surface->GetOutput();
274 }
275
276 vtkSimpleMode::vtkSimpleModeInternal::~vtkSimpleModeInternal()
277 {
278   if(_surface)
279     _surface->Delete();
280 }
281
282 vtkSimpleMode::vtkSimpleMode():Factor(1.),AnimationTime(0.),Internal(new vtkSimpleMode::vtkSimpleModeInternal)
283 {
284   //this->SetInputArrayToProcess(0,0,0,vtkDataObject::FIELD_ASSOCIATION_POINTS,vtkDataSetAttributes::VECTORS);
285 }
286
287 vtkSimpleMode::~vtkSimpleMode()
288 {
289   delete this->Internal;
290 }
291
292 void vtkSimpleMode::SetInputArrayToProcess(int idx, int port, int connection, int ff, const char *name)
293 {
294   if(idx==0)
295     this->Internal->setFieldForReal(name);
296   vtkPolyDataAlgorithm::SetInputArrayToProcess(idx,port,connection,ff,name);
297 }
298
299 double GetOptimalRatioFrom(vtkUnstructuredGrid *dataset, vtkDoubleArray *array)
300 {
301   if(!dataset || !array)
302     throw MZCException("The input dataset and or array is null !");
303   vtkDataArray *coords(dataset->GetPoints()->GetData());
304   vtkDoubleArray *coords2(vtkDoubleArray::SafeDownCast(coords));
305   if(!coords2)
306     throw MZCException("Input coordinates are not float64 !");
307   int nbCompo(array->GetNumberOfComponents());
308   if(coords2->GetNumberOfComponents()!=3 || (nbCompo!=3 && nbCompo!=6))
309     throw MZCException("Input coordinates do not have 3 components as it should !");
310   int nbPts(dataset->GetNumberOfPoints());
311   const double *srcPt1(array->Begin());
312   dataset->ComputeBounds();
313   double *minmax1(dataset->GetBounds());
314   double minmax2[3]={0.,0.,0.};
315   for(int i=0;i<nbPts;i++,srcPt1+=nbCompo)
316     {
317       minmax2[0]=std::max(fabs(srcPt1[0]),minmax2[0]);
318       minmax2[1]=std::max(fabs(srcPt1[1]),minmax2[1]);
319       minmax2[2]=std::max(fabs(srcPt1[2]),minmax2[2]);
320     }
321   double maxDispDelta(*std::max_element(minmax2,minmax2+3));
322   if(maxDispDelta<EPS)
323     maxDispDelta=1.;
324   for(int i=0;i<3;i++)
325     minmax2[i]=minmax1[2*i+1]-minmax1[2*i];
326   double maxGeoDelta(*std::max_element(minmax2,minmax2+3));
327   if(maxDispDelta<EPS)
328     maxDispDelta=1.;
329   return maxGeoDelta/maxDispDelta;
330 }
331
332 int vtkSimpleMode::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
333
334   //std::cerr << "########################################## vtkSimpleMode::RequestInformation ##########################################" << std::endl;
335   try
336     {
337       if(this->Internal->getFieldForReal().empty())
338         return 1;
339       /*vtkUnstructuredGrid *usgIn(0);
340       vtkDoubleArray *arr(0);
341       ExtractInfo(inputVector[0],usgIn,this->Internal->getFieldForReal(),arr);
342       std::vector<std::string> candidatesArrName(GetPossibleArrayNames(usgIn));
343       //
344       double ratio(GetOptimalRatioFrom(usgIn,arr));
345       std::string optArrNameForReal(FindBestRealAmong(candidatesArrName));
346       std::string optArrNameForImag(FindBestImagAmong(candidatesArrName));*/
347       //std::cerr << ratio << std::endl;
348       //std::cerr << optArrNameForReal << " * " << optArrNameForImag << std::endl;
349     }
350   catch(MZCException& e)
351     {
352       std::ostringstream oss;
353       oss << "Exception has been thrown in vtkSimpleMode::RequestInformation : " << e.what() << std::endl;
354       if(this->HasObserver("ErrorEvent") )
355         this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
356       else
357         vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
358       vtkObject::BreakOnError();
359       return 0;
360     }
361   return 1;
362 }
363
364 int vtkSimpleMode::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
365 {
366   //std::cerr << "########################################## vtkSimpleMode::RequestData        ##########################################" << std::endl;
367   try
368     {
369       vtkUnstructuredGrid *usgIn(0);
370       usgIn=ExtractInfo1(inputVector[0]);
371       //
372       int nbPts(usgIn->GetNumberOfPoints());
373       vtkPolyData *outSurface(this->Internal->performConnection(usgIn));
374       vtkDoubleArray *arrRealBase(0);
375       ExtractInfo2(outSurface,this->Internal->getFieldForReal(),arrRealBase);
376       vtkSmartPointer<vtkDoubleArray> arrReal(ForceTo3Compo(arrRealBase));
377       vtkSmartPointer<vtkDoubleArray> arr1(vtkSmartPointer<vtkDoubleArray>::New());
378       arr1->SetName(ZE_DISPLACEMENT_NAME1);
379       arr1->SetNumberOfComponents(3);
380       arr1->SetNumberOfTuples(nbPts);
381       double *ptToFeed1(arr1->Begin());
382       const double *srcPt1(arrReal->Begin());
383       double cst1(Factor*sin(AnimationTime*2*M_PI));
384       std::transform(srcPt1,srcPt1+3*nbPts,ptToFeed1,std::bind2nd(std::multiplies<double>(),cst1));
385       int idx1(outSurface->GetPointData()->AddArray(arr1));
386       outSurface->GetPointData()->SetActiveAttribute(idx1,vtkDataSetAttributes::VECTORS);
387       //
388       //
389       vtkSmartPointer<vtkWarpVector> ws(vtkSmartPointer<vtkWarpVector>::New());//vtkNew
390       ws->SetInputData(outSurface);
391       ws->SetScaleFactor(1.);
392       ws->SetInputArrayToProcess(idx1,0,0,"vtkDataObject::FIELD_ASSOCIATION_POINTS",ZE_DISPLACEMENT_NAME1);
393       ws->Update();
394       vtkSmartPointer<vtkDataSet> ds(ws->GetOutput());
395       ds->GetPointData()->RemoveArray(idx1);
396       //
397       int idx2(0);
398       {
399         vtkDataArray *dummy(0);
400         ExtractInfo3(ds,this->Internal->getFieldForReal(),dummy,idx2);
401       }
402       ds->GetPointData()->SetActiveAttribute(idx2,vtkDataSetAttributes::SCALARS);
403       //
404       vtkInformation *outInfo(outputVector->GetInformationObject(0));
405       vtkPolyData *output(vtkPolyData::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
406       output->ShallowCopy(ds);
407       
408     }
409   catch(MZCException& e)
410     {
411       std::ostringstream oss;
412       oss << "Exception has been thrown in vtkSimpleMode::RequestInformation : " << e.what() << std::endl;
413       if(this->HasObserver("ErrorEvent") )
414         this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
415       else
416         vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
417       vtkObject::BreakOnError();
418       return 0;
419     }
420   return 1;
421 }
422
423 void vtkSimpleMode::PrintSelf(ostream& os, vtkIndent indent)
424 {
425   this->Superclass::PrintSelf(os, indent);
426 }