Salome HOME
47554dc06292d6ada3599f6c5528d7bf640417dd
[modules/visu.git] / src / PIPELINE / VISU_IsoSurfacesPL.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_PipeLine.cxx
24 // Author:  Alexey PETROV
25 // Module : VISU
26 //
27 #include "VISU_IsoSurfacesPL.hxx"
28 #include "VISU_LookupTable.hxx"
29
30 #include "VISU_PipeLineUtils.hxx"
31 #include "VISU_LabelPointsFilter.hxx"
32
33 #include <vtkContourFilter.h>
34
35 #define GAP_COEFFICIENT 0.0001
36
37
38 //----------------------------------------------------------------------------
39 vtkStandardNewMacro(VISU_IsoSurfacesPL);
40
41
42 //----------------------------------------------------------------------------
43 VISU_IsoSurfacesPL
44 ::VISU_IsoSurfacesPL()
45 {
46   SetIsShrinkable(false);
47   SetIsFeatureEdgesAllowed(false);
48
49   SetElnoDisassembleState( true );
50
51   myContourFilter = vtkContourFilter::New();
52
53   myCellDataToPointData = vtkCellDataToPointData::New();
54 }
55
56
57 //----------------------------------------------------------------------------
58 VISU_IsoSurfacesPL
59 ::~VISU_IsoSurfacesPL()
60 {
61   myContourFilter->Delete();
62   myContourFilter = NULL;
63
64   myCellDataToPointData->Delete();
65   myCellDataToPointData = NULL;
66 }
67
68
69 //----------------------------------------------------------------------------
70 unsigned long int 
71 VISU_IsoSurfacesPL
72 ::GetMTime()
73 {
74   unsigned long int aTime = Superclass::GetMTime();
75
76   aTime = std::max(aTime, myCellDataToPointData->GetMTime());
77   aTime = std::max(aTime, myContourFilter->GetMTime());
78
79   return aTime;
80 }
81
82
83 //----------------------------------------------------------------------------
84 void
85 VISU_IsoSurfacesPL
86 ::DoShallowCopy(VISU_PipeLine *thePipeLine,
87                 bool theIsCopyInput)
88 {
89   Superclass::DoShallowCopy(thePipeLine, theIsCopyInput);
90
91   if(VISU_IsoSurfacesPL *aPipeLine = dynamic_cast<VISU_IsoSurfacesPL*>(thePipeLine)){
92     SetNbParts(aPipeLine->GetNbParts());
93     vtkFloatingPointType aRange[2] = {aPipeLine->GetMin(), aPipeLine->GetMax()};
94     SetRange(aRange);
95     SetRangeFixed(aPipeLine->IsRangeFixed());
96   }
97 }
98
99
100 //----------------------------------------------------------------------------
101 int
102 VISU_IsoSurfacesPL
103 ::GetNbParts() 
104 {
105   return myContourFilter->GetNumberOfContours();
106 }
107
108 //----------------------------------------------------------------------------
109 vtkFloatingPointType
110 VISU_IsoSurfacesPL
111 ::GetValue(int i) 
112 {
113   return myContourFilter->GetValue(i);
114 }
115
116
117 //----------------------------------------------------------------------------
118 void
119 VISU_IsoSurfacesPL
120 ::SetNbParts(int theNb) 
121 {
122   myContourFilter->SetNumberOfContours(theNb);
123 }
124
125
126 //----------------------------------------------------------------------------
127 void
128 VISU_IsoSurfacesPL
129 ::SetScalarRange( vtkFloatingPointType theRange[2] ) 
130 {
131   Superclass::SetScalarRange( theRange );
132   SetRange(myRange);
133 }
134
135
136 //----------------------------------------------------------------------------
137 void
138 VISU_IsoSurfacesPL
139 ::SetRange(vtkFloatingPointType theRange[2], bool theIsForced)
140 {
141   if(VISU::CheckIsSameRange(myRange, theRange) && !theIsForced)
142     return;
143
144   if(theRange[0] <= theRange[1]){
145     myRange[0] = theRange[0];  
146     myRange[1] = theRange[1];
147     vtkFloatingPointType aRange[2] = {theRange[0], theRange[1]};
148     if( IsRangeFixed() ) {
149       double aDelta = fabs( aRange[1] - aRange[0] ) * GAP_COEFFICIENT;
150       aRange[0] += aDelta;
151       aRange[1] -= aDelta;
152     }
153     if(GetScaling() == VTK_SCALE_LOG10)
154       VISU_LookupTable::ComputeLogRange(theRange, aRange);
155     myContourFilter->GenerateValues(GetNbParts(), aRange);
156   }
157 }
158
159
160 //----------------------------------------------------------------------------
161 vtkFloatingPointType
162 VISU_IsoSurfacesPL
163 ::GetMin() 
164 {
165   return myRange[0];
166 }
167
168
169 //----------------------------------------------------------------------------
170 vtkFloatingPointType
171 VISU_IsoSurfacesPL
172 ::GetMax() 
173 {
174   return myRange[1];
175 }
176
177
178 //----------------------------------------------------------------------------
179 void
180 VISU_IsoSurfacesPL
181 ::SetRangeFixed(bool theIsFixed)
182 {
183   myIsRangeFixed = theIsFixed;
184   SetRange( myRange, true );
185 }
186
187
188 //----------------------------------------------------------------------------
189 bool
190 VISU_IsoSurfacesPL
191 ::IsRangeFixed()
192 {
193   return myIsRangeFixed;
194 }
195
196
197 void
198 //----------------------------------------------------------------------------
199 VISU_IsoSurfacesPL
200 ::Init()
201 {
202   Superclass::Init();
203
204   SetNbParts(10);
205
206   vtkFloatingPointType aScalarRange[2];
207   GetSourceRange(aScalarRange);
208   SetRange(aScalarRange);
209
210   SetRangeFixed(true);
211 }
212
213 //----------------------------------------------------------------------------
214 void
215 VISU_IsoSurfacesPL
216 ::Build()
217 {
218   Superclass::Build();
219
220   VISU::CellDataToPoint(myContourFilter,
221                         myCellDataToPointData,
222                         GetMergedInput());
223
224 }
225
226
227 //----------------------------------------------------------------------------
228
229 vtkDataSet* 
230 VISU_IsoSurfacesPL
231 ::InsertCustomPL()
232 {
233   return myContourFilter->GetOutput();
234 }
235
236
237 //----------------------------------------------------------------------------
238 unsigned long int
239 VISU_IsoSurfacesPL
240 ::GetMemorySize()
241 {
242   unsigned long int aSize = Superclass::GetMemorySize();
243
244   if(vtkDataSet* aDataSet = myContourFilter->GetOutput())
245     aSize += aDataSet->GetActualMemorySize() * 1024;
246   
247   if(myCellDataToPointData->GetInput())
248     if(vtkDataSet* aDataSet = myCellDataToPointData->GetOutput())
249       aSize += aDataSet->GetActualMemorySize() * 1024;
250
251   return aSize;
252 }
253
254
255 //----------------------------------------------------------------------------
256 void
257 VISU_IsoSurfacesPL
258 ::SetMapScale(vtkFloatingPointType theMapScale)
259 {
260   Superclass::SetMapScale(theMapScale);
261
262   vtkFloatingPointType aRange[2] = {GetMax() - theMapScale*(GetMax()-GetMin()), GetMax()};
263   vtkFloatingPointType aNewRange[2] = {aRange[0], aRange[1]};
264   if( IsRangeFixed() ) {
265     double aDelta = fabs( aNewRange[1] - aNewRange[0] ) * GAP_COEFFICIENT;
266     aNewRange[0] += aDelta;
267     aNewRange[1] -= aDelta;
268   }
269   if(GetScaling() == VTK_SCALE_LOG10)
270     VISU_LookupTable::ComputeLogRange(aRange,aNewRange);
271   myContourFilter->GenerateValues(GetNbParts(), aNewRange);
272 }
273
274
275 //----------------------------------------------------------------------------