Salome HOME
112b55b4a5f1b9c544f08134f24e8672d81e524a
[modules/visu.git] / src / PIPELINE / VISU_DeformedGridPL.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_DeformedGridPL.cxx
24 // Author:  Alexey PETROV
25 // Module : VISU
26 //
27 #include "VISU_DeformedGridPL.hxx"
28 #include "VISU_FieldTransform.hxx"
29 #include "VISU_Plot3DPL.hxx"
30
31 #include "VISU_PipeLineUtils.hxx"
32
33 #include <vtkPolyDataMapper.h>
34 #include <vtkContourFilter.h>
35 #include <vtkWarpScalar.h>
36
37
38 //----------------------------------------------------------------------------
39 vtkStandardNewMacro(VISU_DeformedGridPL);
40
41
42 //----------------------------------------------------------------------------
43 VISU_DeformedGridPL
44 ::VISU_DeformedGridPL():
45   myContourFilter(vtkContourFilter::New()),
46   myWarpScalar(vtkWarpScalar::New()),
47   myIsContour(false),
48   myScaleFactor(1.0),
49   myMapScaleFactor(1.0)
50 {
51   SetIsShrinkable(false);
52   SetNumberOfContours(32);
53 }
54
55
56 //----------------------------------------------------------------------------
57 VISU_DeformedGridPL
58 ::~VISU_DeformedGridPL()
59 {}
60
61
62 //----------------------------------------------------------------------------
63 unsigned long int 
64 VISU_DeformedGridPL
65 ::GetMTime()
66 {
67   unsigned long int aTime = Superclass::GetMTime();
68
69   aTime = std::max(aTime, myContourFilter->GetMTime());
70   aTime = std::max(aTime, myWarpScalar->GetMTime());
71
72   return aTime;
73 }
74
75
76 //----------------------------------------------------------------------------
77 unsigned long int
78 VISU_DeformedGridPL
79 ::GetMemorySize()
80 {
81   unsigned long int aSize = Superclass::GetMemorySize();
82
83   if(vtkDataObject* aDataObject = myContourFilter->GetInput())
84     aSize += aDataObject->GetActualMemorySize() * 1024;
85
86   if(vtkDataObject* aDataObject = myWarpScalar->GetInput())
87     aSize += aDataObject->GetActualMemorySize() * 1024;
88
89   return aSize;
90 }
91
92
93 //----------------------------------------------------------------------------
94 void
95 VISU_DeformedGridPL
96 ::DoShallowCopy(VISU_PipeLine *thePipeLine,
97                 bool theIsCopyInput)
98 {
99   Superclass::DoShallowCopy(thePipeLine, theIsCopyInput);
100
101   if(VISU_DeformedGridPL *aPipeLine = dynamic_cast<VISU_DeformedGridPL*>(thePipeLine)){
102     SetScaleFactor( aPipeLine->GetScaleFactor() );
103     SetContourPrs( aPipeLine->GetIsContourPrs() );
104     SetNumberOfContours( aPipeLine->GetNumberOfContours() );
105   }
106 }
107
108
109 //----------------------------------------------------------------------------
110 void
111 VISU_DeformedGridPL
112 ::Init()
113 {
114   Superclass::Init();
115
116   vtkPointSet* aPointSet = GetFieldTransformFilter()->GetPolyDataOutput();
117   SetScaleFactor( VISU_Plot3DPL::GetScaleFactor( this, aPointSet ) );
118 }
119
120
121 //----------------------------------------------------------------------------
122 void
123 VISU_DeformedGridPL
124 ::Build()
125 {
126   Superclass::Build();
127
128   myWarpScalar->SetInput( GetFieldTransformFilter()->GetPolyDataOutput() );
129   GetPolyDataMapper()->SetInput( myWarpScalar->GetPolyDataOutput() );
130 }
131
132
133 //----------------------------------------------------------------------------
134 void
135 VISU_DeformedGridPL
136 ::Update()
137 {
138   vtkPointSet* aPointSet = GetFieldTransformFilter()->GetPolyDataOutput();
139   if ( !myIsContour ) // surface prs
140   {
141     myWarpScalar->SetInput( aPointSet );
142   }
143   else // contour prs
144   {
145     myContourFilter->SetInput( aPointSet );
146
147     vtkFloatingPointType aScalarRange[2];
148     GetSourceRange( aScalarRange );
149
150     myContourFilter->GenerateValues( GetNumberOfContours(), aScalarRange );
151     myWarpScalar->SetInput( myContourFilter->GetOutput() );
152   }
153
154   Superclass::Update();
155 }
156
157
158 //----------------------------------------------------------------------------
159 void
160 VISU_DeformedGridPL
161 ::SetNumberOfContours(int theNumber)
162 {
163   myContourFilter->SetNumberOfContours(theNumber);
164 }
165
166
167 //----------------------------------------------------------------------------
168 int
169 VISU_DeformedGridPL
170 ::GetNumberOfContours()
171 {
172   return myContourFilter->GetNumberOfContours();
173 }
174
175
176 //----------------------------------------------------------------------------
177 void
178 VISU_DeformedGridPL
179 ::SetScaleFactor(vtkFloatingPointType theScaleFactor)
180 {
181   if ( VISU::CheckIsSameValue( myWarpScalar->GetScaleFactor(), theScaleFactor ) )
182     return;
183
184   myScaleFactor = theScaleFactor;
185   myWarpScalar->SetScaleFactor(theScaleFactor*myMapScaleFactor);
186 }
187
188
189 //----------------------------------------------------------------------------
190 vtkFloatingPointType
191 VISU_DeformedGridPL
192 ::GetScaleFactor()
193 {
194   return myScaleFactor;
195 }
196
197
198 //----------------------------------------------------------------------------
199 void
200 VISU_DeformedGridPL
201 ::SetContourPrs(bool theIsContourPrs )
202 {
203   if(myIsContour == theIsContourPrs)
204     return;
205
206   myIsContour = theIsContourPrs;
207   Modified();
208 }
209
210
211 //----------------------------------------------------------------------------
212 bool
213 VISU_DeformedGridPL
214 ::GetIsContourPrs()
215 {
216   return myIsContour;
217 }
218
219
220 //----------------------------------------------------------------------------
221 void 
222 VISU_DeformedGridPL
223 ::SetMapScale(vtkFloatingPointType theMapScale)
224 {
225   Superclass::SetMapScale(theMapScale);
226   myMapScaleFactor = theMapScale;
227
228   if ( myIsContour ) {
229     vtkFloatingPointType aSourceRange[2];
230     GetSourceRange( aSourceRange );
231     vtkFloatingPointType aDeltaRange = aSourceRange[1] - aSourceRange[0];
232     vtkFloatingPointType aNewRange[2] = { aSourceRange[1] - theMapScale*aDeltaRange, aSourceRange[1] };
233     myContourFilter->GenerateValues( GetNumberOfContours(), aNewRange );
234   }
235
236   myWarpScalar->SetScaleFactor( myScaleFactor * theMapScale );
237 }