Salome HOME
Update copyrights
[modules/paravis.git] / src / Plugins / EllipseBuilder / vtkEllipseBuilderFilter.cxx
1 // Copyright (C) 2014-2019  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 "vtkEllipseBuilderFilter.h"
21
22 #include <vtkObjectFactory.h>
23 #include <vtkInformation.h>
24 #include <vtkInformationVector.h>
25 #include <vtkMultiBlockDataSet.h>
26 #include <vtkDataObjectTreeIterator.h>
27 #include <vtkDataSet.h>
28 #include <vtkPoints.h>
29 #include <vtkPointData.h>
30 #include <vtkDataArray.h>
31 #include <vtkStringArray.h>
32 #include <vtkUnstructuredGrid.h>
33 #include <vtkPolyVertex.h>
34 #include <vtkMath.h>
35
36 #include <complex>
37 #include <algorithm>
38 #include <vector>
39 #include <cmath>
40
41 using namespace std;
42
43 //------------------------------------------------------------------------------
44 bool isStringInList(const list<string>& aList, const string& theName)
45 {
46   list<string>::const_iterator anIter;
47   for (anIter = aList.begin(); anIter != aList.end(); anIter++)
48     {
49       // Compares the values of the string.
50       if ( (*anIter).compare(theName) == 0)
51         return true;
52     }
53   return false;
54 }
55
56 //------------------------------------------------------------------------------
57 void add2List(const list<string>& theSource, list<string>& theDestination)
58 {
59   list<string>::const_iterator anIter;
60   for (anIter = theSource.begin(); anIter != theSource.end(); anIter++)
61     {
62       // Add the item to the list if it does not exist
63       if (!isStringInList(theDestination, *anIter))
64         theDestination.push_back(*anIter);
65     }
66 }
67
68 //------------------------------------------------------------------------------
69 double maximumElement(vector< complex<double> >& theList)
70 {
71   double aMaximum;
72   vector<double> tmpList;
73   vector< complex<double> >::iterator it;
74   for (it = theList.begin(); it != theList.end(); ++it)
75     {
76       tmpList.push_back(((*it)*conj(*it)).real());
77     }
78   aMaximum = *max_element(tmpList.begin(), tmpList.end());
79   return aMaximum;
80 }
81
82 //------------------------------------------------------------------------------
83 list<std::string> GetListOfFields(vtkDataObject* theObject)
84 {
85   list<string> aList;
86   
87   if (theObject->IsA("vtkDataSet"))
88     {
89       vtkDataSet* aDataSet = vtkDataSet::SafeDownCast(theObject);
90       vtkPointData* aPntData = aDataSet->GetPointData();
91       
92       // Add fields name on points
93       int aNbArrays = aPntData->GetNumberOfArrays();
94       for (int i = 0; i<aNbArrays; i++)
95         {
96           if( vtkDataArray* anArray = aPntData->GetArray(i) ) {
97             //Data array should contains at least 3 components (1 - DX, 2 - DY, 3 - DZ)
98             if( anArray->GetNumberOfComponents() >= 3 ) {
99               const char* aName = aPntData->GetArrayName(i);
100               aList.push_back(aName);
101             }
102           }
103         }
104     }
105   return aList;
106 }
107
108 //------------------------------------------------------------------------------
109 vtkStandardNewMacro(vtkEllipseBuilderFilter);
110
111 //------------------------------------------------------------------------------
112 vtkEllipseBuilderFilter::vtkEllipseBuilderFilter() : vtkMultiBlockDataSetAlgorithm()
113 {
114   this->RealField = NULL;
115   this->ImagField = NULL;
116   this->FieldList = vtkStringArray::New();
117 }
118
119 //------------------------------------------------------------------------------
120 vtkEllipseBuilderFilter::~vtkEllipseBuilderFilter()
121 {
122   this->SetRealField(NULL);
123   this->SetRealField(NULL);
124   this->FieldList->Delete();
125 }
126
127 //------------------------------------------------------------------------------
128 int vtkEllipseBuilderFilter::RequestData(vtkInformation* vtkNotUsed(request),
129                                          vtkInformationVector** theInputVector,
130                                          vtkInformationVector* theOutputVector)
131 {
132   int aResolution = this->Resolution;
133   if(aResolution <=0 )
134     return 0;
135   if( this->ScaleFactor == 0.0 )
136     return 0;
137
138   int anAxis = this->Axis;
139   
140   // Get the info objects
141   vtkMultiBlockDataSet* anInputMDSet = vtkMultiBlockDataSet::GetData(theInputVector[0], 0);
142   vtkMultiBlockDataSet* anOutputMDSet = vtkMultiBlockDataSet::GetData(theOutputVector, 0);
143   
144   double a, b, z_point;
145   vector< complex<double> > circle;
146
147   // Points Ellipses
148   double range, min, max, delta;
149   if(this->EndAngle > this->StartAngle) {
150     min = this->StartAngle;
151     max = this->EndAngle;
152   } else {
153     min = this->StartAngle;
154     max = this->EndAngle + 360.;
155   }
156   range = max-min;
157   delta = range/(double)aResolution;
158   for (double iter = min;  iter < max; iter+=delta)
159     {
160       a = cos(vtkMath::RadiansFromDegrees(iter));
161       b = sin(vtkMath::RadiansFromDegrees(iter));
162       complex<double> aVal(a, b);
163       circle.push_back(aVal);
164     }
165   
166   //Scale Factor    
167   double aScaleFactor;
168   aScaleFactor = 1./this->ScaleFactor;  
169     
170   vtkDataObjectTreeIterator* anIter = anInputMDSet->NewTreeIterator();
171   anIter->VisitOnlyLeavesOff();
172   bool created = false;
173   for (anIter->InitTraversal(); !anIter->IsDoneWithTraversal(); anIter->GoToNextItem())
174     {
175       vtkDataObject* anInputNode = anInputMDSet->GetDataSet(anIter);
176       if (anInputNode->IsA("vtkDataSet"))
177         {
178           vtkUnstructuredGrid* aSourceDS = vtkUnstructuredGrid::SafeDownCast(anInputNode);
179           if(!aSourceDS) continue;
180           
181           vtkPointData* aPointData = aSourceDS->GetPointData();
182           if(!aPointData) continue;
183           
184           int aNumberPoints = aSourceDS->GetNumberOfPoints();
185           
186           vtkDataArray* aRealArray = vtkDataArray::SafeDownCast(aPointData->GetArray(this->RealField));
187           vtkDataArray* anImagArray = vtkDataArray::SafeDownCast(aPointData->GetArray(this->ImagField));
188           
189           if(!aRealArray || !anImagArray) continue;
190           
191           int aNumberOfRealComponents = aRealArray->GetNumberOfComponents();
192           int aNumberOfImagComponents = anImagArray->GetNumberOfComponents();              
193           if (aNumberOfRealComponents >= 3 && aNumberOfImagComponents >= 3)
194             {
195               anOutputMDSet->CopyStructure(anInputMDSet);
196               vtkUnstructuredGrid* aCloneDS = aSourceDS->NewInstance();
197               vtkPoints* aClonePoints = vtkPoints::New();
198               aCloneDS->SetPoints(aClonePoints);
199               aClonePoints->Delete();
200               anOutputMDSet->SetDataSet(anIter, aCloneDS);
201               aCloneDS->Delete();             
202               double rx, ry, ix, iy;
203               created = true;
204               for (int j = 0; j < aNumberPoints; j++)
205                 {
206                   z_point = aSourceDS->GetPoint(j)[2];
207                   
208                   if (anAxis == 2)      // Z : DX and DY
209                     {
210                       rx = aRealArray->GetTuple(j)[0];
211                       ry = aRealArray->GetTuple(j)[1];
212                       ix = anImagArray->GetTuple(j)[0];
213                       iy = anImagArray->GetTuple(j)[1];
214                     }
215                   else if (anAxis == 1) // Y : DX and DZ
216                     {
217                       rx = aRealArray->GetTuple(j)[0];
218                       ry = aRealArray->GetTuple(j)[2];
219                       ix = anImagArray->GetTuple(j)[0];
220                       iy = anImagArray->GetTuple(j)[2];
221                     }
222                   else                  // X : DY and DZ
223                     {
224                       rx = aRealArray->GetTuple(j)[1];
225                       ry = aRealArray->GetTuple(j)[2];
226                       ix = anImagArray->GetTuple(j)[1];
227                       iy = anImagArray->GetTuple(j)[2];
228                     }
229                   
230                   complex<double> x(rx, ix);
231                   complex<double> y(ry, iy);               
232                   
233                   x = x / aScaleFactor;
234                   y = y / aScaleFactor;             
235
236                   double x_point, y_point;
237                   for (int r = 0; r < circle.size(); r++)
238                     {
239                       x_point = (x*circle[r]).real();
240                       y_point = (y*circle[r]).real();
241                       vtkIdType anId[1];
242                       if (anAxis == 2)
243                         anId[0] = aClonePoints->InsertNextPoint(x_point, y_point, z_point);
244                       else if (anAxis == 1)
245                         anId[0] = aClonePoints->InsertNextPoint(x_point, z_point, y_point);
246                       else
247                         anId[0] = aClonePoints->InsertNextPoint(z_point, x_point, y_point);                   
248                       aCloneDS->InsertNextCell(VTK_VERTEX, 1, anId);
249                     }             
250                 }
251             }
252           else
253             {
254               continue;
255             }
256         }
257     }
258   anIter->Delete();
259   if(!created)
260     return 0;
261   return 1;
262 }
263
264 //------------------------------------------------------------------------------
265 int vtkEllipseBuilderFilter::RequestInformation(vtkInformation* request,
266                                                 vtkInformationVector **theInputVector,
267                                                 vtkInformationVector *theOutputVector)
268 {
269     // Retrieve an instance of vtkMultiBlockDataSet class from an information object.
270     vtkMultiBlockDataSet* anInputMDataSet = vtkMultiBlockDataSet::GetData(theInputVector[0], 0);
271
272     list<string> aList;
273     vtkDataObjectTreeIterator* anIter = anInputMDataSet->NewTreeIterator();
274     anIter->VisitOnlyLeavesOff();
275     for (anIter->InitTraversal(); !anIter->IsDoneWithTraversal(); anIter->GoToNextItem())
276     {
277         vtkDataObject* anInputNode = anInputMDataSet->GetDataSet(anIter);
278         if (anInputNode->IsA("vtkDataSet"))
279         {
280           list<string> aSubList = GetListOfFields(anInputNode);
281           add2List(aSubList, aList);
282         }
283     }
284     anIter->Delete();
285     
286     this->FieldList->Reset();
287     this->FieldList->SetNumberOfValues(aList.size());
288     list<string>::const_iterator anIterName;
289     int i = 0;
290     for (anIterName = aList.begin(); anIterName != aList.end(); anIterName++)
291     {
292         this->FieldList->SetValue(i, *anIterName);
293         i++;
294     }
295     
296     return this->Superclass::RequestInformation(request, theInputVector, theOutputVector);
297 }
298
299
300 //------------------------------------------------------------------------------
301 vtkStringArray* vtkEllipseBuilderFilter::GetFieldList()
302 {
303     return this->FieldList;
304 }
305
306 //------------------------------------------------------------------------------
307 void vtkEllipseBuilderFilter::PrintSelf(ostream& os, vtkIndent indent)
308 {
309     this->Superclass::PrintSelf(os, indent);
310     os << indent << "Real Field : " << this->RealField << endl;
311     os << indent << "Imag Field : " << this->ImagField << endl;
312     os << indent << "Start Angle : " << this->StartAngle << endl;
313     os << indent << "End Angle : " << this->EndAngle << endl;
314     os << indent << "Scale Factor : " << this->ScaleFactor << endl;
315     os << indent << "Scale Resolution : " << this->Resolution << endl;
316     os << indent << "Axis : " << this->Axis << endl;
317 }
318