Salome HOME
Fix for Bug IPAL8945
[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::VISU_ScalarMapOnDeformedShapePL(){
52   myDeformVectors = vtkWarpVector::New();
53   myMergeFilter   = vtkMergeFilter::New();
54   myExtractorScalars = VISU_Extractor::New();
55   myCellDataToPointData = vtkCellDataToPointData::New();
56 }
57
58 /*!
59  * Destructor.
60  * Delete all fields.
61 */
62 VISU_ScalarMapOnDeformedShapePL::~VISU_ScalarMapOnDeformedShapePL(){
63   myDeformVectors->UnRegisterAllOutputs();
64   myDeformVectors->Delete();
65
66   myMergeFilter->UnRegisterAllOutputs();
67   myMergeFilter->Delete();
68   
69   myExtractorScalars->UnRegisterAllOutputs();
70   myExtractorScalars->Delete();
71
72   myCellDataToPointData->UnRegisterAllOutputs();
73   myCellDataToPointData->Delete();
74 }
75
76 /*!
77  * Initial method
78  */
79 void VISU_ScalarMapOnDeformedShapePL::Init(){
80
81   if (GetScalars() == NULL) SetScalars(GetInput2());
82   
83   Superclass::Init();
84   
85   float aScalarRange[2];
86   GetSourceRange(aScalarRange);
87   static double EPS = 1.0 / VTK_LARGE_FLOAT;
88   if(aScalarRange[1] > EPS)
89     SetScale(VISU_DeformedShapePL::GetScaleFactor(GetInput2())/aScalarRange[1]);
90   else
91     SetScale(0.0);
92
93   myMapper->SetColorModeToMapScalars();
94   myMapper->ScalarVisibilityOn();
95
96   // Sets input for field transformation filter
97   myFieldTransform->SetInput(myExtractor->GetOutput());
98
99 }
100
101 /*!
102  * Build method
103  * Building of deformation and puts result to merge filter.
104  */
105 void VISU_ScalarMapOnDeformedShapePL::Build()
106 {
107   // Set input for extractor
108   myExtractor->SetInput(GetInput2());
109   
110   VISU::CellDataToPoint(myDeformVectors,myCellDataToPointData,
111                         GetInput2(),myFieldTransform);
112   
113   // Sets geometry for merge filter
114   myMergeFilter->SetGeometry(myDeformVectors->GetOutput());
115   // Sets data to mapper
116   myMapper->SetInput(myMergeFilter->GetOutput());
117 }
118
119 /*!
120  *  Update method
121  */
122 void VISU_ScalarMapOnDeformedShapePL::Update(){
123   this->UpdateScalars();
124   
125   float aRange[2];
126   GetSourceRange(aRange);
127   float aScalarRange[2] = {aRange[0], aRange[1]};
128   
129   if(myBarTable->GetScale() == VTK_SCALE_LOG10)
130     VISU_LookupTable::ComputeLogRange(aRange,aScalarRange);
131   myMapperTable->SetRange(aScalarRange);
132
133   myMapperTable->Build();
134   myBarTable->Build();
135
136   myMapper->SetLookupTable(myMapperTable);
137   myMapper->SetScalarRange(aScalarRange);
138   
139   VISU_PipeLine::Update();
140 }
141
142 /*!
143  * Update scalars method.
144  * Put scalars to merge filter.
145  */
146 void VISU_ScalarMapOnDeformedShapePL::UpdateScalars(){
147   if(myScalars->GetCellData()->GetVectors() != NULL ||
148      myScalars->GetPointData()->GetVectors() != NULL)
149     myMergeFilter->SetScalars(myExtractorScalars->GetOutput());
150   else
151     myMergeFilter->SetScalars(GetScalars());
152 }
153
154 /*!
155  * Copy information about pipline.
156  * Copy scale and scalars.
157  */
158 void VISU_ScalarMapOnDeformedShapePL::ShallowCopy(VISU_PipeLine *thePipeLine){
159   VISU_ScalarMapOnDeformedShapePL *aPipeLine = dynamic_cast<VISU_ScalarMapOnDeformedShapePL*>(thePipeLine);
160   if(aPipeLine){
161      SetScale(aPipeLine->GetScale());
162      SetScalars(aPipeLine->GetScalars());
163      float aRange[2];
164      aPipeLine->GetSourceRange(aRange);
165      SetScalarRange(aRange);
166   }
167   Superclass::ShallowCopy(thePipeLine);
168 }
169
170 /*!
171  * Set scalars.
172  * Sets vtkDataSet with scalars values to VISU_Extractor filter for scalars extraction.
173  */
174 void VISU_ScalarMapOnDeformedShapePL::SetScalars(vtkDataSet *theScalars){
175   myScalars = theScalars;
176   vtkUnstructuredGrid* aScalars = GetScalars();
177   myExtractorScalars->SetInput(aScalars);
178   Modified();
179 }
180
181 /*!
182  * Get pointer to input scalars.
183  */
184 vtkUnstructuredGrid* VISU_ScalarMapOnDeformedShapePL::GetScalars(){
185   return myScalars.GetPointer();
186 }
187
188 /*!
189  * Sets scale for deformed shape
190  */
191 void VISU_ScalarMapOnDeformedShapePL::SetScale(float theScale) {
192   if(myScaleFactor == theScale) return;
193   myScaleFactor = theScale;
194   myDeformVectors->SetScaleFactor(myScaleFactor);
195   Modified();
196 }
197
198 /*!
199  * Gets scale of deformed shape.
200  */
201 float VISU_ScalarMapOnDeformedShapePL::GetScale() {
202   float aScale=myDeformVectors->GetScaleFactor();
203   return aScale;
204 }
205
206 /*!
207  * Set scale factor of deformation.
208  */
209 void VISU_ScalarMapOnDeformedShapePL::SetMapScale(float theMapScale){
210   myDeformVectors->SetScaleFactor(myScaleFactor*theMapScale);
211   Modified();
212 }
213
214 /*!
215  * Gets scalar mode.
216  */
217 int VISU_ScalarMapOnDeformedShapePL::GetScalarMode(){
218   int aMode=myExtractorScalars->GetScalarMode();
219   return aMode;
220 }
221
222 /*!
223  * Sets scalar mode.
224  */
225 void VISU_ScalarMapOnDeformedShapePL::SetScalarMode(int theScalarMode){
226   myExtractorScalars->SetScalarMode(theScalarMode);
227   Modified();
228 }
229
230 /*!
231  * Gets ranges of extracted scalars
232  * \param theRange[2] - output values
233  * \li theRange[0] - minimum value
234  * \li theRange[1] - maximum value
235  */
236 void VISU_ScalarMapOnDeformedShapePL::GetSourceRange(float theRange[2]){
237   myExtractorScalars->Update();
238   myExtractorScalars->GetUnstructuredGridOutput()->GetScalarRange(theRange);
239 }