Salome HOME
IMP 0017328: min and max scalar map of results given at gauss points. A fix for the...
[modules/visu.git] / src / PIPELINE / VISU_FieldTransform.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 //  File   : VISU_FieldTransform.cxx
23 //  Module : VISU
24 //
25 #include "VISU_FieldTransform.hxx"
26 #include "VTKViewer_Transform.h"
27 #include "VISU_PipeLineUtils.hxx"
28
29 #include <vtkObjectFactory.h>
30 #include <vtkPointData.h>
31 #include <vtkCellData.h>
32 #include <vtkDataSet.h>
33 #include <vtkMath.h>
34 #include <vtkInformation.h>
35 #include <vtkInformationVector.h>
36
37 static vtkFloatingPointType Tolerance = 1.0 / VTK_LARGE_FLOAT;
38
39 //----------------------------------------------------------------------------
40 vtkStandardNewMacro(VISU_FieldTransform);
41
42 //----------------------------------------------------------------------------
43 double
44 VISU_FieldTransform
45 ::Ident(double theArg)
46 {
47   return theArg;
48 }
49
50 //----------------------------------------------------------------------------
51 double
52 VISU_FieldTransform
53 ::Log10(double theArg)
54 {
55   if(theArg <= 0.0) 
56     return -VTK_LARGE_FLOAT;
57
58   return log10(theArg);
59 }
60
61
62 //----------------------------------------------------------------------------
63 VISU_FieldTransform
64 ::VISU_FieldTransform()
65 {
66   myFunction = &Ident;
67   myTransform = NULL;
68
69   myScalarRange[0] = VTK_LARGE_FLOAT;
70   myScalarRange[1] = -VTK_LARGE_FLOAT;
71 }
72
73 //----------------------------------------------------------------------------
74 VISU_FieldTransform
75 ::~VISU_FieldTransform() 
76 {
77   SetSpaceTransform(NULL);
78 }
79
80
81 unsigned long 
82 VISU_FieldTransform
83 ::GetMTime()
84 {
85   unsigned long aTime = Superclass::GetMTime();
86   if(myTransform)
87     aTime = std::max(aTime, myTransform->GetMTime());
88
89   return aTime;
90 }
91
92 //----------------------------------------------------------------------------
93 void
94 VISU_FieldTransform
95 ::SetScalarTransform(TTransformFun theFunction) 
96 {
97   if(myFunction == theFunction)
98     return;
99
100   if(theFunction == NULL) 
101     theFunction = &Ident;
102
103   myFunction = theFunction;
104
105   Modified();
106 }
107
108 //----------------------------------------------------------------------------
109 void 
110 VISU_FieldTransform
111 ::SetSpaceTransform(VTKViewer_Transform* theTransform)
112 {
113   if(myTransform == theTransform)
114     return;
115
116   if(myTransform != NULL) 
117     myTransform->UnRegister(this);
118
119   myTransform = theTransform;
120
121   if(theTransform != NULL) 
122     theTransform->Register(this);
123
124   Modified();
125 }
126
127
128 //----------------------------------------------------------------------------
129 void
130 VISU_FieldTransform
131 ::SetScalarRange(vtkFloatingPointType theScalarRange[2]) 
132 {
133   if(VISU::CheckIsSameRange(theScalarRange, myScalarRange))
134     return;
135
136   myScalarRange[0] = theScalarRange[0];
137   myScalarRange[1] = theScalarRange[1];
138
139   Modified();
140 }
141
142 //----------------------------------------------------------------------------
143 void
144 VISU_FieldTransform
145 ::SetScalarMin(vtkFloatingPointType theValue)
146 {
147   vtkFloatingPointType aScalarRange[2] = {theValue, GetScalarRange()[1]};
148   SetScalarRange(aScalarRange);
149 }
150
151 //----------------------------------------------------------------------------
152 void
153 VISU_FieldTransform
154 ::SetScalarMax(vtkFloatingPointType theValue)
155 {
156   vtkFloatingPointType aScalarRange[2] = {GetScalarRange()[0], theValue};
157   SetScalarRange(aScalarRange);
158 }
159
160 //----------------------------------------------------------------------------
161 template<typename TValueType> 
162 void
163 LinearTransformVectors(TValueType* theInputPtr,
164                        TValueType* theOutputPtr,
165                        vtkIdType theNbOfTuples,
166                        vtkFloatingPointType theScale[3])
167 {
168   for(vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++){
169     for(vtkIdType aComponentId = 0; aComponentId < 3; aComponentId++){
170       *theOutputPtr = TValueType(*theInputPtr * theScale[aComponentId]);
171       theOutputPtr++;
172       theInputPtr++;
173     }
174   }
175 }
176
177
178 //----------------------------------------------------------------------------
179 template<typename TValueType> 
180 void
181 NonLinearTransformVectors(vtkDataArray *theInputVectors,
182                           TValueType* theInputPtr,
183                           TValueType* theOutputPtr,
184                           vtkIdType theNbOfTuples,
185                           vtkFloatingPointType theScale[3],
186                           VISU_FieldTransform::TTransformFun theFunction,
187                           vtkFloatingPointType theModifiedScalarMin,
188                           vtkFloatingPointType theModifiedScalarDelta,
189                           vtkFloatingPointType theSourceScalarMax)
190 {
191   for(vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++){
192     vtkFloatingPointType anInputVector[3];
193     theInputVectors->GetTuple(aTupleId, anInputVector);
194     vtkFloatingPointType aMagnification = vtkMath::Norm(anInputVector);
195     if(aMagnification > Tolerance)
196       aMagnification = 
197         ((*theFunction)(aMagnification) - theModifiedScalarMin) / 
198         theModifiedScalarDelta * theSourceScalarMax / 
199         aMagnification;
200     if(aMagnification < 0.0) 
201       aMagnification = 0.0;
202     for(vtkIdType aComponentId = 0; aComponentId < 3; aComponentId++){
203       *theOutputPtr = TValueType(*theInputPtr * aMagnification * theScale[aComponentId]);
204       theOutputPtr++;
205       theInputPtr++;
206     }
207   }
208 }
209
210
211 //----------------------------------------------------------------------------
212 template<typename TDataSetAttributesType> 
213 void
214 ExecuteVectors(VISU_FieldTransform::TTransformFun theFunction,
215                VTKViewer_Transform* theTransform,
216                vtkFloatingPointType theScalarRange[2], 
217                vtkIdType theNbOfTuples,
218                TDataSetAttributesType* theInputData, 
219                TDataSetAttributesType* theOutputData)
220 {
221   vtkDataArray *anInputVectors = theInputData->GetVectors();
222   if(!anInputVectors || theNbOfTuples < 1) 
223     return;
224
225   vtkFloatingPointType aScalarRange[2];
226   aScalarRange[0] = (*theFunction)(theScalarRange[0]);
227   aScalarRange[1] = (*theFunction)(theScalarRange[1]);
228
229   vtkFloatingPointType aScalarDelta = aScalarRange[1] - aScalarRange[0];
230   vtkFloatingPointType aScale[3] = {1.0, 1.0, 1.0};
231
232   if(theTransform){
233     aScale[0] = theTransform->GetScale()[0];
234     aScale[1] = theTransform->GetScale()[1];
235     aScale[2] = theTransform->GetScale()[2];
236   }
237
238   vtkIdType anInputDataType = anInputVectors->GetDataType();
239   vtkDataArray *anOutputVectors = vtkDataArray::CreateDataArray(anInputDataType);
240   anOutputVectors->SetNumberOfComponents(3);
241   anOutputVectors->SetNumberOfTuples(theNbOfTuples);
242
243   void *anInputPtr = anInputVectors->GetVoidPointer(0);
244   void *anOutputPtr = anOutputVectors->GetVoidPointer(0);
245
246   if(theFunction == &(VISU_FieldTransform::Ident)){
247     switch(anInputDataType){
248       vtkTemplateMacro4(LinearTransformVectors,
249                         (VTK_TT *)(anInputPtr), 
250                         (VTK_TT *)(anOutputPtr), 
251                         theNbOfTuples,
252                         aScale);
253     default:
254       break;
255     }  
256   }else{
257     switch(anInputDataType){
258       vtkTemplateMacro9(NonLinearTransformVectors,
259                         anInputVectors,
260                         (VTK_TT *)(anInputPtr), 
261                         (VTK_TT *)(anOutputPtr), 
262                         theNbOfTuples,
263                         aScale,
264                         theFunction,
265                         aScalarRange[0],
266                         aScalarDelta,
267                         theScalarRange[1]);
268     default:
269       break;
270     }  
271   }
272
273   theOutputData->SetVectors(anOutputVectors);
274   anOutputVectors->Delete();
275 }
276
277
278 //----------------------------------------------------------------------------
279 template<typename TValueType> 
280 void
281 NonLinearTransformScalars(vtkDataArray *theInputScalars,
282                           TValueType* theInputPtr,
283                           TValueType* theOutputPtr,
284                           vtkIdType theNbOfTuples,
285                           VISU_FieldTransform::TTransformFun theFunction,
286                           vtkFloatingPointType theModifiedScalarMin)
287 {
288   for(vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++){
289     vtkFloatingPointType aScalar = (*theFunction)(vtkFloatingPointType(*theInputPtr));
290     if(aScalar < theModifiedScalarMin) 
291       aScalar = theModifiedScalarMin;
292     *theOutputPtr = TValueType(aScalar);
293     theOutputPtr++;
294     theInputPtr++;
295   }
296 }
297
298
299 //----------------------------------------------------------------------------
300 template<typename TDataSetAttributesType> 
301 void
302 ExecuteScalars(VISU_FieldTransform::TTransformFun theFunction, 
303                vtkFloatingPointType theScalarRange[2],
304                vtkIdType theNbOfTuples, 
305                TDataSetAttributesType* theInputData, 
306                TDataSetAttributesType* theOutputData)
307 {
308   vtkDataArray *anInputScalars = theInputData->GetScalars();
309   if(!anInputScalars || theNbOfTuples < 1)
310     return;
311
312   vtkFloatingPointType aScalarRange[2];
313   aScalarRange[0] = (*theFunction)(theScalarRange[0]);
314   aScalarRange[1] = (*theFunction)(theScalarRange[1]);
315
316   vtkIdType anInputDataType = anInputScalars->GetDataType();
317   vtkDataArray *anOutputScalars = vtkDataArray::CreateDataArray(anInputDataType);
318   anOutputScalars->SetNumberOfComponents(1);
319   anOutputScalars->SetNumberOfTuples(theNbOfTuples);
320
321   void *anInputPtr = anInputScalars->GetVoidPointer(0);
322   void *anOutputPtr = anOutputScalars->GetVoidPointer(0);
323
324   switch(anInputDataType){
325     vtkTemplateMacro6(NonLinearTransformScalars,
326                       anInputScalars,
327                       (VTK_TT *)(anInputPtr), 
328                       (VTK_TT *)(anOutputPtr), 
329                       theNbOfTuples,
330                       theFunction,
331                       aScalarRange[0]);
332   default:
333     break;
334   }  
335   
336   theOutputData->SetScalars(anOutputScalars);
337   anOutputScalars->Delete();
338 }
339
340
341 //----------------------------------------------------------------------------
342 int
343 VISU_FieldTransform
344 ::RequestData(vtkInformation *vtkNotUsed(request),
345               vtkInformationVector **inputVector,
346               vtkInformationVector *outputVector)
347 {
348   // get the info objects
349   vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
350   vtkInformation *outInfo = outputVector->GetInformationObject(0);
351
352   // get the input and ouptut
353   vtkDataSet *input = vtkDataSet::SafeDownCast(
354     inInfo->Get(vtkDataObject::DATA_OBJECT()));
355   vtkDataSet *output = vtkDataSet::SafeDownCast(
356     outInfo->Get(vtkDataObject::DATA_OBJECT()));
357
358   output->CopyStructure(input);
359   if(myFunction != &Ident || (myTransform && !myTransform->IsIdentity())){
360     output->GetPointData()->CopyScalarsOff();
361     output->GetPointData()->CopyVectorsOff();
362
363     output->GetCellData()->CopyScalarsOff();
364     output->GetCellData()->CopyVectorsOff();
365
366     ExecuteScalars(myFunction,
367                    myScalarRange,
368                    input->GetNumberOfPoints(),
369                    input->GetPointData(),
370                    output->GetPointData());
371
372     ExecuteVectors(myFunction,
373                    myTransform,
374                    myScalarRange,
375                    input->GetNumberOfPoints(),
376                    input->GetPointData(),
377                    output->GetPointData());
378
379     ExecuteScalars(myFunction,
380                    myScalarRange,
381                    input->GetNumberOfCells(),
382                    input->GetCellData(),
383                    output->GetCellData());
384
385     ExecuteVectors(myFunction,
386                    myTransform,
387                    myScalarRange,
388                    input->GetNumberOfCells(),
389                    input->GetCellData(),
390                    output->GetCellData());
391   }else{
392     output->GetPointData()->CopyAllOn();
393     output->GetCellData()->CopyAllOn();
394
395     output->GetPointData()->PassData(input->GetPointData());
396     output->GetCellData()->PassData(input->GetCellData());
397   }
398
399   output->GetPointData()->PassData(input->GetPointData());
400   output->GetCellData()->PassData(input->GetCellData());
401   return 1;
402 }
403
404
405 //----------------------------------------------------------------------------