Salome HOME
Update copyright information
[modules/visu.git] / src / PIPELINE / VISU_Extractor.cxx
1 // Copyright (C) 2007-2012  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
10 //
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 // Lesser General Public License for more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22
23 //  VISU OBJECT : interactive object for VISU entities implementation
24 //  File   : VISU_Extractor.cxx
25 //  Module : VISU
26 //
27 #include "VISU_Extractor.hxx"
28 #include "VISU_PipeLineUtils.hxx"
29 #include "VISU_ConvertorUtils.hxx"
30 #include "VISU_MeshValue.hxx"
31
32 #include <SUIT_Session.h>
33 #include <SUIT_ResourceMgr.h>
34
35 #include <sstream>
36 #include <cmath> // to use std::abs( int )
37
38 #include <vtkObjectFactory.h>
39 #include <vtkUnstructuredGrid.h>
40 #include <vtkPointData.h>
41 #include <vtkCellData.h>
42 #include <vector>
43 #include <vtkInformation.h>
44 #include <vtkInformationVector.h>
45
46 // Algorithm to cut values (see CutValue())
47 // Use sprintf function
48 #define USE_SPRINTF
49 // Use new algorithm (works only when USE_SPRINTF is disabled)
50 #define USE_NEW_ALGORITHM
51
52 //----------------------------------------------------------------------------
53 vtkStandardNewMacro(VISU_Extractor);
54
55 //----------------------------------------------------------------------------
56 VISU_Extractor
57 ::VISU_Extractor():
58   myScalarMode(1),
59   myGaussMetric(VISU::AVERAGE_METRIC)
60 {
61 }
62
63
64 //----------------------------------------------------------------------------
65 VISU_Extractor
66 ::~VISU_Extractor()
67 {}
68
69
70 //----------------------------------------------------------------------------
71 void
72 VISU_Extractor
73 ::SetScalarMode(int theScalarMode)
74 {
75   if(myScalarMode != theScalarMode){
76     myScalarMode = theScalarMode;
77     Modified();
78   }
79 }
80
81 //----------------------------------------------------------------------------
82 int
83 VISU_Extractor
84 ::GetScalarMode()
85 {
86   return myScalarMode;
87 }
88
89
90 //----------------------------------------------------------------------------
91 void
92 VISU_Extractor
93 ::SetGaussMetric(VISU::TGaussMetric theGaussMetric)
94 {
95   if(myGaussMetric != theGaussMetric){
96     myGaussMetric = theGaussMetric;
97     Modified();
98   }
99 }
100
101 //----------------------------------------------------------------------------
102 VISU::TGaussMetric
103 VISU_Extractor
104 ::GetGaussMetric()
105 {
106   return myGaussMetric;
107 }
108
109 //----------------------------------------------------------------------------
110 vtkFloatingPointType CutValue (vtkFloatingPointType theValue, int theDecimals)
111 {
112   vtkFloatingPointType v = theValue;
113
114 #if defined(USE_SPRINTF)
115   char aFormat[16];
116   sprintf(aFormat, "%%.%dg", theDecimals);
117   char aStr [256];
118   sprintf(aStr, aFormat, theValue);
119   v = atof(aStr);
120 #elif defined(USE_NEW_ALGORITHM)
121   //
122   // VSR 19/10/2009: new algorithm does not use long long type
123   //
124   long n1 = 0;
125   // calculate order of the integral part
126   while ( v > 1. ) {
127     v = v / 10.;
128     n1++;
129   }
130   // calculate length of the fractional part
131   long n2 = theDecimals - n1;
132   v = theValue;
133   if ( n2 > 0 )
134     while ( n2-- > 0 ) v = v * 10.;
135   else
136     while ( n2++ < 0 ) v = v / 10.;
137   v = floor( (double) v );
138   n2 = theDecimals - n1;
139   if ( n2 > 0 )
140     while ( n2-- > 0 ) v = v / 10.;
141   else
142     while ( n2++ < 0 ) v = v * 10.;
143 #else
144   //
145   // VSR 19/10/2009: old algorithm uses long long type -
146   // causes incompatibility with some platforms (Windows?)
147   //
148   vtkFloatingPointType aDegree = 0.0;
149   if (abs((long long)v) > 1)
150     aDegree = (long long)log10((double)abs((long long)v)) + 1;
151   aDegree = theDecimals - aDegree;
152   //printf("$$$ 1 v = %.20g , aDegree = %lld \n", v, (long long)aDegree);
153
154   aDegree = pow(10, aDegree);
155   v = ((vtkFloatingPointType)((long long)(v * aDegree))) / aDegree;
156   //printf("$$$ 2 v = %.20g , aDegree = %lld \n", v, (long long)aDegree);
157 #endif
158
159   return v;
160 }
161
162 //----------------------------------------------------------------------------
163 template<typename TValueType>
164 void
165 Module2Scalars(vtkDataArray *theInputDataArray,
166                TValueType* theOutputPtr,
167                vtkIdType theNbOfTuples)
168 {
169   vtkIdType aNbComp = theInputDataArray->GetNumberOfComponents();
170   std::vector<vtkFloatingPointType> anArray(aNbComp < 3? 3: aNbComp);
171   for(vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++){
172     theInputDataArray->GetTuple(aTupleId, &anArray[0]);
173     vtkFloatingPointType aVector[3] = {anArray[0], anArray[1], anArray[2]};
174     vtkFloatingPointType aScalar = sqrt(aVector[0]*aVector[0] +
175                                         aVector[1]*aVector[1] +
176                                         aVector[2]*aVector[2]);
177     *theOutputPtr = TValueType(aScalar);
178     theOutputPtr++;
179   }
180 }
181
182 template<typename TValueType>
183 void
184 Module2ScalarsMOD(vtkDataArray *theInputDataArray,
185                   TValueType* theOutputPtr,
186                   vtkIdType theNbOfTuples,
187                   VISU::TGaussMetric theGaussMetric)
188 {
189   vtkIdType aNbComp = theInputDataArray->GetNumberOfComponents();
190   if (aNbComp != 3) // Min, Max, Avg
191     return;
192   std::vector<vtkFloatingPointType> anArray (3);
193   for (vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++) {
194     theInputDataArray->GetTuple(aTupleId, &anArray[0]);
195     switch (theGaussMetric) {
196       case VISU::MINIMUM_METRIC: *theOutputPtr = TValueType(anArray[0]); break;
197       case VISU::MAXIMUM_METRIC: *theOutputPtr = TValueType(anArray[1]); break;
198       case VISU::AVERAGE_METRIC: *theOutputPtr = TValueType(anArray[2]); break;
199     }
200     theOutputPtr++;
201   }
202 }
203
204
205 //----------------------------------------------------------------------------
206 template<typename TValueType>
207 void
208 Component2Scalars(vtkDataArray *theInputDataArray,
209                   TValueType* theInputPtr,
210                   TValueType* theOutputPtr,
211                   vtkIdType theNbOfTuples,
212                   vtkIdType theComponentId)
213 {
214   vtkIdType aNbComp = theInputDataArray->GetNumberOfComponents();
215   for (vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++) {
216     *theOutputPtr = *(theInputPtr + theComponentId);
217     theInputPtr += aNbComp;
218     theOutputPtr++;
219   }
220 }
221
222 //----------------------------------------------------------------------------
223 template<typename TValueType>
224 void
225 CutScalarsTempl(TValueType* theDataPtr,
226                 vtkIdType theNbOfTuples,
227                 int theDecimals)
228 {
229   for (vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++) {
230     *theDataPtr = TValueType(CutValue(*theDataPtr, theDecimals));
231     theDataPtr++;
232   }
233 }
234
235 //----------------------------------------------------------------------------
236 template<typename TDataSetAttributesType> void
237 ExecuteScalars(vtkIdType theNbOfTuples,
238                vtkIdType theScalarMode,
239                VISU::TGaussMetric theGaussMetric,
240                TDataSetAttributesType* theInputData,
241                TDataSetAttributesType* theOutputData)
242 {
243   if (theNbOfTuples < 1)
244     return;
245
246   vtkDataArray* aFieldArray = NULL;
247   switch (theGaussMetric) {
248     case VISU::AVERAGE_METRIC: aFieldArray = theInputData->GetArray("VISU_FIELD"); break;
249     case VISU::MINIMUM_METRIC: aFieldArray = theInputData->GetArray("VISU_FIELD_GAUSS_MIN"); break;
250     case VISU::MAXIMUM_METRIC: aFieldArray = theInputData->GetArray("VISU_FIELD_GAUSS_MAX"); break;
251   }
252   if( !aFieldArray )
253     return;
254
255   vtkIdType anInputDataType = aFieldArray->GetDataType();
256   vtkDataArray *anOutputScalars = vtkDataArray::CreateDataArray(anInputDataType);
257   anOutputScalars->SetNumberOfComponents(1);
258   anOutputScalars->SetNumberOfTuples(theNbOfTuples);
259
260   void *anInputPtr = aFieldArray->GetVoidPointer(0);
261   void *anOutputPtr = anOutputScalars->GetVoidPointer(0);
262
263   if (theScalarMode == 0) {
264     vtkDataArray* aFieldArrayMOD = theInputData->GetArray("VISU_FIELD_GAUSS_MOD");
265     if (aFieldArrayMOD) {
266       switch (anInputDataType) {
267         vtkTemplateMacro4(Module2ScalarsMOD,
268                           aFieldArrayMOD,
269                           (VTK_TT *)(anOutputPtr),
270                           theNbOfTuples,
271                           theGaussMetric);
272       default:
273         break;
274       }
275     }
276     else {
277       switch (anInputDataType) {
278         vtkTemplateMacro3(Module2Scalars,
279                           aFieldArray,
280                           (VTK_TT *)(anOutputPtr),
281                           theNbOfTuples);
282       default:
283         break;
284       }
285     }
286   } else {
287     switch (anInputDataType) {
288       vtkTemplateMacro5(Component2Scalars,
289                         aFieldArray,
290                         (VTK_TT *)(anInputPtr),
291                         (VTK_TT *)(anOutputPtr),
292                         theNbOfTuples,
293                         theScalarMode - 1);
294     default:
295       break;
296     }
297   }
298
299   theOutputData->SetScalars(anOutputScalars);
300   anOutputScalars->Delete();
301 }
302
303 //---------------------------------------------------------------
304 template<typename TDataSetAttributesType> void
305 CutScalars(vtkIdType theNbOfTuples,
306            TDataSetAttributesType* theData)
307 {
308   if (theNbOfTuples < 1)
309     return;
310
311   vtkDataArray *aScalars = theData->GetScalars();
312   if (!aScalars)
313     return;
314
315   vtkIdType aDataType = aScalars->GetDataType();
316   void *aPtr = aScalars->GetVoidPointer(0);
317
318   SUIT_ResourceMgr* aResourceMgr = SUIT_Session::session()->resourceMgr();
319   // san: precision can be negative - this is used by double spin boxes
320   int aDecimals = std::abs( aResourceMgr->integerValue("VISU", "visual_data_precision", 6) );
321
322   switch(aDataType) {
323     vtkTemplateMacro3(CutScalarsTempl,
324                       (VTK_TT *)(aPtr),
325                       theNbOfTuples,
326                       aDecimals);
327   default:
328     break;
329   }
330 }
331
332 //---------------------------------------------------------------
333 int
334 VISU_Extractor
335 ::RequestData(vtkInformation *theRequest,
336               vtkInformationVector **theInputVector,
337               vtkInformationVector *theOutputVector)
338 {
339   vtkDataSet *anInput = VISU::GetInput( theInputVector, 0 );
340   vtkDataSet *anOutput = VISU::GetOutput( theOutputVector );
341
342   anOutput->CopyStructure( anInput );
343
344   vtkPointData *anInputPointData = anInput->GetPointData();
345   vtkPointData *anOutputPointData = anOutput->GetPointData();
346   anOutputPointData->PassData( anInputPointData );
347   if ( VISU::IsDataOnPoints( anInput ) ) {
348     int aNbElems = anInput->GetNumberOfPoints();
349     if ( anInputPointData->GetAttribute( vtkDataSetAttributes::VECTORS ) )
350       ExecuteScalars( aNbElems, myScalarMode, myGaussMetric, anInputPointData, anOutputPointData );
351     CutScalars( aNbElems, anOutputPointData );
352   }
353
354   vtkCellData *anInputCellData = anInput->GetCellData();
355   vtkCellData *anOutputCellData = anOutput->GetCellData();
356   anOutputCellData->PassData( anInputCellData );
357   if ( VISU::IsDataOnCells( anInput ) ) {
358     int aNbElems = anInput->GetNumberOfCells();
359     if ( anInputCellData->GetAttribute( vtkDataSetAttributes::VECTORS ) )
360       ExecuteScalars( aNbElems, myScalarMode, myGaussMetric, anInputCellData, anOutputCellData );
361     CutScalars( aNbElems, anOutputCellData );
362   }
363
364   return 1;
365 }