]> SALOME platform Git repositories - modules/visu.git/blob - src/PIPELINE/VISU_ScalarMapOnDeformedShapePL.cxx
Salome HOME
Join modifications from branch OCC_debug_for_3_2_0b1
[modules/visu.git] / src / PIPELINE / VISU_ScalarMapOnDeformedShapePL.cxx
1 //  VISU ScalarMapOnDeformedShapePL
2 //
3 //  Copyright (C) 2003  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.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org
21 //
22 //
23 // File:    VISU_ScalarMapOnDeformedShapePL.cxx
24 // Author:  Eugeny Nikolaev
25 // Module : VISU
26
27 #include "VISU_ScalarMapOnDeformedShapePL.hxx"
28 #include "VISU_DeformedShapePL.hxx"
29 #include "VISU_PipeLineUtils.hxx"
30 #include "VTKViewer_TransformFilter.h"
31 #include "VTKViewer_Transform.h"
32
33 #include <vtkWarpVector.h>
34 #include <vtkMergeFilter.h>
35 #include <vtkUnstructuredGrid.h>
36 #include <vtkCellDataToPointData.h>
37 #include <vtkPointDataToCellData.h>
38
39 #include <VISU_Convertor.hxx>
40
41 vtkStandardNewMacro(VISU_ScalarMapOnDeformedShapePL)
42
43 /*!
44  * Constructor. Creating new instances of vtkWarpVector,vtkMergeFilter,vtkUnstructuredGrid
45  * Where:
46  * \li myDeformVectors is vtkWarpVector  - deformation vectors data
47  * \li myMergeFilter   is vtkMergeFilter - merge filter.
48  * Merge filter which unify the deformation and scalars
49  * \li myScalars is vtk shared pointer to vtkUnstructuredGrid - scalars data
50 */
51 VISU_ScalarMapOnDeformedShapePL
52 ::VISU_ScalarMapOnDeformedShapePL()
53 {
54   myDeformVectors = vtkWarpVector::New();
55   myMergeFilter   = vtkMergeFilter::New();
56   myExtractorScalars = VISU_Extractor::New();
57   myCellDataToPointData = vtkCellDataToPointData::New();
58 }
59
60 /*!
61  * Destructor.
62  * Delete all fields.
63 */
64 VISU_ScalarMapOnDeformedShapePL
65 ::~VISU_ScalarMapOnDeformedShapePL()
66 {
67   myDeformVectors->UnRegisterAllOutputs();
68   myDeformVectors->Delete();
69
70   myMergeFilter->UnRegisterAllOutputs();
71   myMergeFilter->Delete();
72   
73   myExtractorScalars->UnRegisterAllOutputs();
74   myExtractorScalars->Delete();
75
76   myCellDataToPointData->UnRegisterAllOutputs();
77   myCellDataToPointData->Delete();
78 }
79
80 /*!
81  * Initial method
82  */
83 void
84 VISU_ScalarMapOnDeformedShapePL
85 ::Init()
86 {
87   if (GetScalars() == NULL) SetScalars(GetInput2());
88   
89   Superclass::Init();
90   
91   vtkFloatingPointType aScalarRange[2];
92   GetSourceRange(aScalarRange);
93   static double EPS = 1.0 / VTK_LARGE_FLOAT;
94   if(aScalarRange[1] > EPS)
95     SetScale(VISU_DeformedShapePL::GetScaleFactor(GetInput2())/aScalarRange[1]);
96   else
97     SetScale(0.0);
98
99   myMapper->SetColorModeToMapScalars();
100   myMapper->ScalarVisibilityOn();
101
102   // Sets input for field transformation filter
103   myFieldTransform->SetInput(myExtractor->GetOutput());
104
105 }
106
107 /*!
108  * Build method
109  * Building of deformation and puts result to merge filter.
110  */
111 void
112 VISU_ScalarMapOnDeformedShapePL
113 ::Build()
114 {
115   // Set input for extractor
116   myExtractor->SetInput(GetInput2());
117   
118   VISU::CellDataToPoint(myDeformVectors,myCellDataToPointData,
119                         GetInput2(),myFieldTransform);
120   
121   // Sets geometry for merge filter
122   myMergeFilter->SetGeometry(myDeformVectors->GetOutput());
123   // Sets data to mapper
124   myMapper->SetInput(myMergeFilter->GetOutput());
125 }
126
127 /*!
128  *  Update method
129  */
130 void
131 VISU_ScalarMapOnDeformedShapePL
132 ::Update()
133 {
134   this->UpdateScalars();
135   
136   vtkFloatingPointType aRange[2];
137   GetSourceRange(aRange);
138   vtkFloatingPointType aScalarRange[2] = {aRange[0], aRange[1]};
139   
140   if(myBarTable->GetScale() == VTK_SCALE_LOG10)
141     VISU_LookupTable::ComputeLogRange(aRange,aScalarRange);
142   myMapperTable->SetRange(aScalarRange);
143
144   myMapperTable->Build();
145   myBarTable->Build();
146
147   myMapper->SetLookupTable(myMapperTable);
148   myMapper->SetScalarRange(aScalarRange);
149   
150   VISU_PipeLine::Update();
151 }
152
153 /*!
154  * Update scalars method.
155  * Put scalars to merge filter.
156  */
157 void
158 VISU_ScalarMapOnDeformedShapePL
159 ::UpdateScalars()
160 {
161   if(myScalars->GetCellData()->GetVectors() != NULL ||
162      myScalars->GetPointData()->GetVectors() != NULL)
163     myMergeFilter->SetScalars(myExtractorScalars->GetOutput());
164   else
165     myMergeFilter->SetScalars(GetScalars());
166 }
167
168 /*!
169  * Copy information about pipline.
170  * Copy scale and scalars.
171  */
172 void
173 VISU_ScalarMapOnDeformedShapePL
174 ::ShallowCopy(VISU_PipeLine *thePipeLine)
175 {
176   VISU_ScalarMapOnDeformedShapePL *aPipeLine = dynamic_cast<VISU_ScalarMapOnDeformedShapePL*>(thePipeLine);
177   if(aPipeLine){
178      SetScale(aPipeLine->GetScale());
179      SetScalars(aPipeLine->GetScalars());
180      vtkFloatingPointType aRange[2];
181      aPipeLine->GetSourceRange(aRange);
182      SetScalarRange(aRange);
183   }
184   Superclass::ShallowCopy(thePipeLine);
185 }
186
187 /*!
188  * Set scalars.
189  * Sets vtkDataSet with scalars values to VISU_Extractor filter for scalars extraction.
190  */
191 void
192 VISU_ScalarMapOnDeformedShapePL
193 ::SetScalars(vtkDataSet *theScalars)
194 {
195   myScalars = theScalars;
196   vtkUnstructuredGrid* aScalars = GetScalars();
197   myExtractorScalars->SetInput(aScalars);
198   Modified();
199 }
200
201 /*!
202  * Get pointer to input scalars.
203  */
204 vtkUnstructuredGrid* 
205 VISU_ScalarMapOnDeformedShapePL
206 ::GetScalars()
207 {
208   return myScalars.GetPointer();
209 }
210
211 /*!
212  * Sets scale for deformed shape
213  */
214 void
215 VISU_ScalarMapOnDeformedShapePL
216 ::SetScale(vtkFloatingPointType theScale) 
217 {
218   if(myScaleFactor == theScale) return;
219   myScaleFactor = theScale;
220   myDeformVectors->SetScaleFactor(myScaleFactor);
221   Modified();
222 }
223
224 /*!
225  * Gets scale of deformed shape.
226  */
227 vtkFloatingPointType
228 VISU_ScalarMapOnDeformedShapePL
229 ::GetScale() 
230 {
231   vtkFloatingPointType aScale=myDeformVectors->GetScaleFactor();
232   return aScale;
233 }
234
235 /*!
236  * Set scale factor of deformation.
237  */
238 void
239 VISU_ScalarMapOnDeformedShapePL
240 ::SetMapScale(vtkFloatingPointType theMapScale)
241 {
242   myDeformVectors->SetScaleFactor(myScaleFactor*theMapScale);
243   Modified();
244 }
245
246 /*!
247  * Gets scalar mode.
248  */
249 int
250 VISU_ScalarMapOnDeformedShapePL
251 ::GetScalarMode()
252 {
253   int aMode=myExtractorScalars->GetScalarMode();
254   return aMode;
255 }
256
257 /*!
258  * Sets scalar mode.
259  */
260 void
261 VISU_ScalarMapOnDeformedShapePL
262 ::SetScalarMode(int theScalarMode)
263 {
264   myExtractorScalars->SetScalarMode(theScalarMode);
265   Modified();
266 }
267
268 /*!
269  * Gets ranges of extracted scalars
270  * \param theRange[2] - output values
271  * \li theRange[0] - minimum value
272  * \li theRange[1] - maximum value
273  */
274 void 
275 VISU_ScalarMapOnDeformedShapePL
276 ::GetSourceRange(vtkFloatingPointType theRange[2])
277 {
278   myExtractorScalars->Update();
279   myExtractorScalars->GetUnstructuredGridOutput()->GetScalarRange(theRange);
280 }