Salome HOME
1dbbeb0596eb0f2d12854bfa187d7d92bd106403
[modules/visu.git] / src / PIPELINE / VISU_DeformedGridPL.cxx
1 // Copyright (C) 2007-2011  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 //  VISU OBJECT : interactive object for VISU entities implementation
21 // File:    VISU_DeformedGridPL.cxx
22 // Author:  Alexey PETROV
23 // Module : VISU
24 //
25 #include "VISU_DeformedGridPL.hxx"
26 #include "VISU_FieldTransform.hxx"
27 #include "VISU_Plot3DPL.hxx"
28
29 #include "VISU_PipeLineUtils.hxx"
30
31 #include <vtkPolyDataMapper.h>
32 #include <vtkContourFilter.h>
33 #include <vtkWarpScalar.h>
34
35
36 //----------------------------------------------------------------------------
37 vtkStandardNewMacro(VISU_DeformedGridPL);
38
39
40 //----------------------------------------------------------------------------
41 VISU_DeformedGridPL
42 ::VISU_DeformedGridPL():
43   myContourFilter(vtkContourFilter::New()),
44   myWarpScalar(vtkWarpScalar::New()),
45   myIsContour(false),
46   myScaleFactor(1.0),
47   myMapScaleFactor(1.0)
48 {
49   SetIsShrinkable(false);
50   SetNumberOfContours(32);
51 }
52
53
54 //----------------------------------------------------------------------------
55 VISU_DeformedGridPL
56 ::~VISU_DeformedGridPL()
57 {}
58
59
60 //----------------------------------------------------------------------------
61 unsigned long int 
62 VISU_DeformedGridPL
63 ::GetMTime()
64 {
65   unsigned long int aTime = Superclass::GetMTime();
66
67   aTime = std::max(aTime, myContourFilter->GetMTime());
68   aTime = std::max(aTime, myWarpScalar->GetMTime());
69
70   return aTime;
71 }
72
73
74 //----------------------------------------------------------------------------
75 unsigned long int
76 VISU_DeformedGridPL
77 ::GetMemorySize()
78 {
79   unsigned long int aSize = Superclass::GetMemorySize();
80
81   if(vtkDataObject* aDataObject = myContourFilter->GetInput())
82     aSize += aDataObject->GetActualMemorySize() * 1024;
83
84   if(vtkDataObject* aDataObject = myWarpScalar->GetInput())
85     aSize += aDataObject->GetActualMemorySize() * 1024;
86
87   return aSize;
88 }
89
90
91 //----------------------------------------------------------------------------
92 void
93 VISU_DeformedGridPL
94 ::DoShallowCopy(VISU_PipeLine *thePipeLine,
95                 bool theIsCopyInput)
96 {
97   Superclass::DoShallowCopy(thePipeLine, theIsCopyInput);
98
99   if(VISU_DeformedGridPL *aPipeLine = dynamic_cast<VISU_DeformedGridPL*>(thePipeLine)){
100     SetScaleFactor( aPipeLine->GetScaleFactor() );
101     SetContourPrs( aPipeLine->GetIsContourPrs() );
102     SetNumberOfContours( aPipeLine->GetNumberOfContours() );
103   }
104 }
105
106
107 //----------------------------------------------------------------------------
108 void
109 VISU_DeformedGridPL
110 ::Init()
111 {
112   Superclass::Init();
113
114   vtkPointSet* aPointSet = GetFieldTransformFilter()->GetPolyDataOutput();
115   SetScaleFactor( VISU_Plot3DPL::GetScaleFactor( this, aPointSet ) );
116 }
117
118
119 //----------------------------------------------------------------------------
120 void
121 VISU_DeformedGridPL
122 ::Build()
123 {
124   Superclass::Build();
125
126   myWarpScalar->SetInput( GetFieldTransformFilter()->GetPolyDataOutput() );
127   GetPolyDataMapper()->SetInput( myWarpScalar->GetPolyDataOutput() );
128 }
129
130
131 //----------------------------------------------------------------------------
132 void
133 VISU_DeformedGridPL
134 ::Update()
135 {
136   vtkPointSet* aPointSet = GetFieldTransformFilter()->GetPolyDataOutput();
137   if ( !myIsContour ) // surface prs
138   {
139     myWarpScalar->SetInput( aPointSet );
140   }
141   else // contour prs
142   {
143     myContourFilter->SetInput( aPointSet );
144
145     vtkFloatingPointType aScalarRange[2];
146     GetSourceRange( aScalarRange );
147
148     myContourFilter->GenerateValues( GetNumberOfContours(), aScalarRange );
149     myWarpScalar->SetInput( myContourFilter->GetOutput() );
150   }
151
152   Superclass::Update();
153 }
154
155
156 //----------------------------------------------------------------------------
157 void
158 VISU_DeformedGridPL
159 ::SetNumberOfContours(int theNumber)
160 {
161   myContourFilter->SetNumberOfContours(theNumber);
162 }
163
164
165 //----------------------------------------------------------------------------
166 int
167 VISU_DeformedGridPL
168 ::GetNumberOfContours()
169 {
170   return myContourFilter->GetNumberOfContours();
171 }
172
173
174 //----------------------------------------------------------------------------
175 void
176 VISU_DeformedGridPL
177 ::SetScaleFactor(vtkFloatingPointType theScaleFactor)
178 {
179   if ( VISU::CheckIsSameValue( myWarpScalar->GetScaleFactor(), theScaleFactor ) )
180     return;
181
182   myScaleFactor = theScaleFactor;
183   myWarpScalar->SetScaleFactor(theScaleFactor*myMapScaleFactor);
184 }
185
186
187 //----------------------------------------------------------------------------
188 vtkFloatingPointType
189 VISU_DeformedGridPL
190 ::GetScaleFactor()
191 {
192   return myScaleFactor;
193 }
194
195
196 //----------------------------------------------------------------------------
197 void
198 VISU_DeformedGridPL
199 ::SetContourPrs(bool theIsContourPrs )
200 {
201   if(myIsContour == theIsContourPrs)
202     return;
203
204   myIsContour = theIsContourPrs;
205   Modified();
206 }
207
208
209 //----------------------------------------------------------------------------
210 bool
211 VISU_DeformedGridPL
212 ::GetIsContourPrs()
213 {
214   return myIsContour;
215 }
216
217
218 //----------------------------------------------------------------------------
219 void 
220 VISU_DeformedGridPL
221 ::SetMapScale(vtkFloatingPointType theMapScale)
222 {
223   Superclass::SetMapScale(theMapScale);
224   myMapScaleFactor = theMapScale;
225
226   if ( myIsContour ) {
227     vtkFloatingPointType aSourceRange[2];
228     GetSourceRange( aSourceRange );
229     vtkFloatingPointType aDeltaRange = aSourceRange[1] - aSourceRange[0];
230     vtkFloatingPointType aNewRange[2] = { aSourceRange[1] - theMapScale*aDeltaRange, aSourceRange[1] };
231     myContourFilter->GenerateValues( GetNumberOfContours(), aNewRange );
232   }
233
234   myWarpScalar->SetScaleFactor( myScaleFactor * theMapScale );
235 }