]> SALOME platform Git repositories - modules/visu.git/blob - src/PIPELINE/VISU_Extractor.cxx
Salome HOME
52831311df4204ac3eaeaab31a830990be8e73bc
[modules/visu.git] / src / PIPELINE / VISU_Extractor.cxx
1 //  Copyright (C) 2007-2008  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 //  VISU OBJECT : interactive object for VISU entities implementation
23 //  File   : VISU_Extractor.cxx
24 //  Module : VISU
25
26 #include "VISU_Extractor.hxx"
27 #include "VISU_PipeLineUtils.hxx"
28 #include "VISU_ConvertorUtils.hxx"
29 #include "VISU_MeshValue.hxx"
30
31 #include <SUIT_Session.h>
32 #include <SUIT_ResourceMgr.h>
33
34 #include <sstream>
35
36 #include <vtkObjectFactory.h>
37 #include <vtkUnstructuredGrid.h>
38 #include <vtkPointData.h>
39 #include <vtkCellData.h>
40 #include <vector>
41 #include <vtkInformation.h>
42 #include <vtkInformationVector.h>
43
44 #define USE_SPRINTF 1
45
46 //----------------------------------------------------------------------------
47 vtkStandardNewMacro(VISU_Extractor);
48
49 //----------------------------------------------------------------------------
50 VISU_Extractor
51 ::VISU_Extractor():
52   myScalarMode(1),
53   myGaussMetric(VISU::AVERAGE_METRIC)
54 {
55 }
56
57
58 //----------------------------------------------------------------------------
59 VISU_Extractor
60 ::~VISU_Extractor()
61 {}
62
63
64 //----------------------------------------------------------------------------
65 void
66 VISU_Extractor
67 ::SetScalarMode(int theScalarMode)
68 {
69   if(myScalarMode != theScalarMode){
70     myScalarMode = theScalarMode;
71     Modified();
72   }
73 }
74
75 //----------------------------------------------------------------------------
76 int
77 VISU_Extractor
78 ::GetScalarMode()
79 {
80   return myScalarMode;
81 }
82
83
84 //----------------------------------------------------------------------------
85 void
86 VISU_Extractor
87 ::SetGaussMetric(VISU::TGaussMetric theGaussMetric)
88 {
89   if(myGaussMetric != theGaussMetric){
90     myGaussMetric = theGaussMetric;
91     Modified();
92   }
93 }
94
95 //----------------------------------------------------------------------------
96 VISU::TGaussMetric
97 VISU_Extractor
98 ::GetGaussMetric()
99 {
100   return myGaussMetric;
101 }
102
103 //----------------------------------------------------------------------------
104 vtkFloatingPointType CutValue (vtkFloatingPointType theValue, int theDecimals)
105 {
106   vtkFloatingPointType v = theValue;
107
108 #ifdef USE_SPRINTF
109   char aFormat[16];
110   sprintf(aFormat, "%%.%dg", theDecimals);
111   char aStr [256];
112   sprintf(aStr, aFormat, theValue);
113   v = atof(aStr);
114 #else
115
116 #ifndef USE_OLD_ALGORITHM
117   //
118   // VSR 19/10/2009: new algorithm does not use long long type
119   //
120   long n1 = 0;
121   // calculate order of the integral part
122   while ( v > 1. ) {
123     v = v / 10.;
124     n1++;
125   }
126   // calculate length of the fractional part
127   long n2 = theDecimals - n1;
128   v = theValue;
129   if ( n2 > 0 )
130     while ( n2-- > 0 ) v = v * 10.;
131   else
132     while ( n2++ < 0 ) v = v / 10.;
133   v = floor( (double) v );
134   n2 = theDecimals - n1;
135   if ( n2 > 0 )
136     while ( n2-- > 0 ) v = v / 10.;
137   else
138     while ( n2++ < 0 ) v = v * 10.;
139 #else
140   //
141   // VSR 19/10/2009: old algorithm uses long long type -
142   // causes incompatibility with some platforms (Windows?)
143   //
144   vtkFloatingPointType aDegree = 0.0;
145   if (abs((long long)v) > 1)
146     aDegree = (long long)log10((double)abs((long long)v)) + 1;
147   aDegree = theDecimals - aDegree;
148   //printf("$$$ 1 v = %.20g , aDegree = %lld \n", v, (long long)aDegree);
149
150   aDegree = pow(10, aDegree);
151   v = ((vtkFloatingPointType)((long long)(v * aDegree))) / aDegree;
152   //printf("$$$ 2 v = %.20g , aDegree = %lld \n", v, (long long)aDegree);
153 #endif
154
155 #endif
156
157   return v;
158 }
159
160 //----------------------------------------------------------------------------
161 template<typename TValueType>
162 void
163 Module2Scalars(vtkDataArray *theInputDataArray,
164                TValueType* theOutputPtr,
165                vtkIdType theNbOfTuples)
166 {
167   vtkIdType aNbComp = theInputDataArray->GetNumberOfComponents();
168   std::vector<vtkFloatingPointType> anArray(aNbComp < 3? 3: aNbComp);
169   for(vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++){
170     theInputDataArray->GetTuple(aTupleId, &anArray[0]);
171     vtkFloatingPointType aVector[3] = {anArray[0], anArray[1], anArray[2]};
172     vtkFloatingPointType aScalar = sqrt(aVector[0]*aVector[0] +
173                                         aVector[1]*aVector[1] +
174                                         aVector[2]*aVector[2]);
175     *theOutputPtr = TValueType(aScalar);
176     theOutputPtr++;
177   }
178 }
179
180 template<typename TValueType>
181 void
182 Module2ScalarsMOD(vtkDataArray *theInputDataArray,
183                   TValueType* theOutputPtr,
184                   vtkIdType theNbOfTuples,
185                   VISU::TGaussMetric theGaussMetric)
186 {
187   vtkIdType aNbComp = theInputDataArray->GetNumberOfComponents();
188   if (aNbComp != 3) // Min, Max, Avg
189     return;
190   std::vector<vtkFloatingPointType> anArray (3);
191   for (vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++) {
192     theInputDataArray->GetTuple(aTupleId, &anArray[0]);
193     switch (theGaussMetric) {
194       case VISU::MINIMUM_METRIC: *theOutputPtr = TValueType(anArray[0]); break;
195       case VISU::MAXIMUM_METRIC: *theOutputPtr = TValueType(anArray[1]); break;
196       case VISU::AVERAGE_METRIC: *theOutputPtr = TValueType(anArray[2]); break;
197     }
198     theOutputPtr++;
199   }
200 }
201
202
203 //----------------------------------------------------------------------------
204 template<typename TValueType>
205 void
206 Component2Scalars(vtkDataArray *theInputDataArray,
207                   TValueType* theInputPtr,
208                   TValueType* theOutputPtr,
209                   vtkIdType theNbOfTuples,
210                   vtkIdType theComponentId)
211 {
212   vtkIdType aNbComp = theInputDataArray->GetNumberOfComponents();
213   for (vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++) {
214     *theOutputPtr = *(theInputPtr + theComponentId);
215     theInputPtr += aNbComp;
216     theOutputPtr++;
217   }
218 }
219
220 //----------------------------------------------------------------------------
221 template<typename TValueType>
222 void
223 CutScalarsTempl(TValueType* theDataPtr,
224                 vtkIdType theNbOfTuples,
225                 int theDecimals)
226 {
227   for (vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++) {
228     *theDataPtr = TValueType(CutValue(*theDataPtr, theDecimals));
229     theDataPtr++;
230   }
231 }
232
233 //----------------------------------------------------------------------------
234 template<typename TDataSetAttributesType> void
235 ExecuteScalars(vtkIdType theNbOfTuples,
236                vtkIdType theScalarMode,
237                VISU::TGaussMetric theGaussMetric,
238                TDataSetAttributesType* theInputData,
239                TDataSetAttributesType* theOutputData)
240 {
241   if (theNbOfTuples < 1)
242     return;
243
244   vtkDataArray* aFieldArray = NULL;
245   switch (theGaussMetric) {
246     case VISU::AVERAGE_METRIC: aFieldArray = theInputData->GetArray("VISU_FIELD"); break;
247     case VISU::MINIMUM_METRIC: aFieldArray = theInputData->GetArray("VISU_FIELD_GAUSS_MIN"); break;
248     case VISU::MAXIMUM_METRIC: aFieldArray = theInputData->GetArray("VISU_FIELD_GAUSS_MAX"); break;
249   }
250   if( !aFieldArray )
251     return;
252
253   vtkIdType anInputDataType = aFieldArray->GetDataType();
254   vtkDataArray *anOutputScalars = vtkDataArray::CreateDataArray(anInputDataType);
255   anOutputScalars->SetNumberOfComponents(1);
256   anOutputScalars->SetNumberOfTuples(theNbOfTuples);
257
258   void *anInputPtr = aFieldArray->GetVoidPointer(0);
259   void *anOutputPtr = anOutputScalars->GetVoidPointer(0);
260
261   if (theScalarMode == 0) {
262     vtkDataArray* aFieldArrayMOD = theInputData->GetArray("VISU_FIELD_GAUSS_MOD");
263     if (aFieldArrayMOD) {
264       switch (anInputDataType) {
265         vtkTemplateMacro4(Module2ScalarsMOD,
266                           aFieldArrayMOD,
267                           (VTK_TT *)(anOutputPtr),
268                           theNbOfTuples,
269                           theGaussMetric);
270       default:
271         break;
272       }
273     }
274     else {
275       switch (anInputDataType) {
276         vtkTemplateMacro3(Module2Scalars,
277                           aFieldArray,
278                           (VTK_TT *)(anOutputPtr),
279                           theNbOfTuples);
280       default:
281         break;
282       }
283     }
284   } else {
285     switch (anInputDataType) {
286       vtkTemplateMacro5(Component2Scalars,
287                         aFieldArray,
288                         (VTK_TT *)(anInputPtr),
289                         (VTK_TT *)(anOutputPtr),
290                         theNbOfTuples,
291                         theScalarMode - 1);
292     default:
293       break;
294     }
295   }
296
297   theOutputData->SetScalars(anOutputScalars);
298   anOutputScalars->Delete();
299 }
300
301 //---------------------------------------------------------------
302 template<typename TDataSetAttributesType> void
303 CutScalars(vtkIdType theNbOfTuples,
304            TDataSetAttributesType* theData)
305 {
306   if (theNbOfTuples < 1)
307     return;
308
309   vtkDataArray *aScalars = theData->GetScalars();
310   if (!aScalars)
311     return;
312
313   vtkIdType aDataType = aScalars->GetDataType();
314   void *aPtr = aScalars->GetVoidPointer(0);
315
316   SUIT_ResourceMgr* aResourceMgr = SUIT_Session::session()->resourceMgr();
317   int aDecimals = aResourceMgr->integerValue("VISU", "floating_point_precision", 6);
318
319   switch(aDataType) {
320     vtkTemplateMacro3(CutScalarsTempl,
321                       (VTK_TT *)(aPtr),
322                       theNbOfTuples,
323                       aDecimals);
324   default:
325     break;
326   }
327 }
328
329 //---------------------------------------------------------------
330 int
331 VISU_Extractor
332 ::RequestData(vtkInformation *theRequest,
333               vtkInformationVector **theInputVector,
334               vtkInformationVector *theOutputVector)
335 {
336   vtkDataSet *anInput = VISU::GetInput( theInputVector, 0 );
337   vtkDataSet *anOutput = VISU::GetOutput( theOutputVector );
338
339   anOutput->CopyStructure( anInput );
340
341   vtkPointData *anInputPointData = anInput->GetPointData();
342   vtkPointData *anOutputPointData = anOutput->GetPointData();
343   anOutputPointData->PassData( anInputPointData );
344   if ( VISU::IsDataOnPoints( anInput ) ) {
345     int aNbElems = anInput->GetNumberOfPoints();
346     if ( anInputPointData->GetAttribute( vtkDataSetAttributes::VECTORS ) )
347       ExecuteScalars( aNbElems, myScalarMode, myGaussMetric, anInputPointData, anOutputPointData );
348     CutScalars( aNbElems, anOutputPointData );
349   }
350
351   vtkCellData *anInputCellData = anInput->GetCellData();
352   vtkCellData *anOutputCellData = anOutput->GetCellData();
353   anOutputCellData->PassData( anInputCellData );
354   if ( VISU::IsDataOnCells( anInput ) ) {
355     int aNbElems = anInput->GetNumberOfCells();
356     if ( anInputCellData->GetAttribute( vtkDataSetAttributes::VECTORS ) )
357       ExecuteScalars( aNbElems, myScalarMode, myGaussMetric, anInputCellData, anOutputCellData );
358     CutScalars( aNbElems, anOutputCellData );
359   }
360
361   return 1;
362 }