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