Salome HOME
790b38a59a5a6d194a46c5b5d7f7f71e8ec123a0
[modules/visu.git] / src / PIPELINE / VISU_VectorsPL.cxx
1 // Copyright (C) 2007-2011  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
23 //  VISU OBJECT : interactive object for VISU entities implementation
24 // File:    VISU_PipeLine.cxx
25 // Author:  Alexey PETROV
26 // Module : VISU
27 //
28 #include "VISU_VectorsPL.hxx"
29 #include "VISU_FieldTransform.hxx"
30 #include "VISU_PipeLineUtils.hxx"
31 #include "VISU_UsedPointsFilter.hxx"
32 #include "VTKViewer_TransformFilter.h"
33 #include "VTKViewer_Transform.h"
34 #include "VTKViewer_CellCenters.h"
35
36 #include <vtkGlyph3D.h>
37 #include <vtkConeSource.h>
38 #include <vtkLineSource.h>
39 #include <vtkGlyphSource2D.h>
40 #include <vtkPolyData.h>
41
42
43 //----------------------------------------------------------------------------
44 vtkStandardNewMacro(VISU_VectorsPL);
45
46 //----------------------------------------------------------------------------
47 template<class TOutputFilter>
48 void ToCellCenters( TOutputFilter *theOutputFilter, 
49                     VTKViewer_CellCenters *theCellCenters,
50                     vtkDataSet* theDataSet,
51                     VISU_UsedPointsFilter* theUsedPointsFilter )
52 {
53   if ( VISU::IsDataOnCells( theDataSet ) ) {
54     theCellCenters->SetInput( theDataSet );
55     theCellCenters->VertexCellsOn();
56     theOutputFilter->SetInput( theCellCenters->GetOutput() );
57   }else {
58     theUsedPointsFilter->SetInput( theDataSet );
59     theOutputFilter->SetInput( theUsedPointsFilter->GetOutput() );
60   }
61 }
62
63
64 VISU_VectorsPL
65 ::VISU_VectorsPL()
66 {
67   SetIsShrinkable(false);
68   SetIsFeatureEdgesAllowed(false);
69
70   myGlyph = vtkGlyph3D::New();
71
72   myGlyphSource = vtkGlyphSource2D::New();
73   myConeSource = vtkConeSource::New();
74   myLineSource = vtkLineSource::New();
75
76   myCenters = VTKViewer_CellCenters::New();
77   myTransformFilter = VTKViewer_TransformFilter::New();
78
79   myUsedPointsFilter = VISU_UsedPointsFilter::New();
80 }
81
82
83 //----------------------------------------------------------------------------
84 VISU_VectorsPL
85 ::~VISU_VectorsPL()
86 {
87   myGlyph->Delete();
88
89   myCenters->Delete();
90
91   myGlyphSource->Delete();
92
93   myConeSource->Delete();
94
95   myLineSource->Delete();
96
97   myTransformFilter->Delete();
98   
99   myUsedPointsFilter->Delete();
100 }
101
102
103 //----------------------------------------------------------------------------
104 unsigned long int 
105 VISU_VectorsPL
106 ::GetMTime()
107 {
108   unsigned long int aTime = Superclass::GetMTime();
109
110   aTime = std::max(aTime, myGlyph->GetMTime());
111   aTime = std::max(aTime, myCenters->GetMTime());
112   aTime = std::max(aTime, myGlyphSource->GetMTime());
113   aTime = std::max(aTime, myConeSource->GetMTime());
114   aTime = std::max(aTime, myLineSource->GetMTime());
115   aTime = std::max(aTime, myTransformFilter->GetMTime());
116
117   return aTime;
118 }
119
120
121 //----------------------------------------------------------------------------
122 void
123 VISU_VectorsPL
124 ::DoShallowCopy(VISU_PipeLine *thePipeLine,
125                 bool theIsCopyInput)
126 {
127   Superclass::DoShallowCopy(thePipeLine, theIsCopyInput);
128
129   if(VISU_VectorsPL *aPipeLine = dynamic_cast<VISU_VectorsPL*>(thePipeLine)){
130     SetGlyphType(aPipeLine->GetGlyphType());
131     SetGlyphPos(aPipeLine->GetGlyphPos());
132   }
133 }
134
135
136 //----------------------------------------------------------------------------
137 void
138 VISU_VectorsPL
139 ::SetTransform( VTKViewer_Transform* theTransform )
140 {
141   myTransformFilter->SetTransform( theTransform );
142   myTransformFilter->Modified();
143 }
144
145
146 //----------------------------------------------------------------------------
147 VTKViewer_Transform* 
148 VISU_VectorsPL
149 ::GetTransform()
150 {
151   return static_cast< VTKViewer_Transform* >( myTransformFilter->GetTransform() );
152 }
153
154
155 //----------------------------------------------------------------------------
156 void
157 VISU_VectorsPL
158 ::SetScale(vtkFloatingPointType theScale) 
159 {
160   if ( VISU::CheckIsSameValue( myScaleFactor, theScale ) )
161     return;
162
163   myGlyph->SetScaleFactor( theScale );
164
165   myScaleFactor = theScale;
166
167   Modified();
168 }
169
170
171 //----------------------------------------------------------------------------
172 vtkFloatingPointType
173 VISU_VectorsPL
174 ::GetScale() 
175 {
176   return myGlyph->GetScaleFactor();
177 }
178
179
180 //----------------------------------------------------------------------------
181 void
182 VISU_VectorsPL
183 ::SetGlyphType(VISU_VectorsPL::GlyphType theType) 
184 {
185   if(myTypeGlyph == theType)
186     return;
187
188   myTypeGlyph = theType;
189   Modified();
190 }
191
192
193 //----------------------------------------------------------------------------
194 VISU_VectorsPL::GlyphType
195 VISU_VectorsPL
196 ::GetGlyphType() const
197 {
198   return myTypeGlyph;
199 }
200
201
202 //----------------------------------------------------------------------------
203 void
204 VISU_VectorsPL
205 ::SetGlyphPos(VISU_VectorsPL::GlyphPos thePos) 
206 {
207   if(myPosGlyph == thePos)
208     return;
209
210   myPosGlyph = thePos;
211   Modified();
212 }
213
214
215 //----------------------------------------------------------------------------
216 VISU_VectorsPL::GlyphPos
217 VISU_VectorsPL
218 ::GetGlyphPos() const
219 {
220   return myPosGlyph;
221 }
222
223
224 //----------------------------------------------------------------------------
225 void
226 VISU_VectorsPL
227 ::Init()
228 {
229   Superclass::Init();
230
231   SetGlyphType(ARROW);
232   SetGlyphPos(TAIL);
233 }
234
235
236 //----------------------------------------------------------------------------
237 void
238 VISU_VectorsPL
239 ::Build()
240 {
241   Superclass::Build();
242   
243   ToCellCenters( myTransformFilter,
244                  myCenters,
245                  GetMergedInput(),
246                  myUsedPointsFilter );
247
248   myGlyph->SetInput( myTransformFilter->GetOutput() );
249   myGlyph->SetVectorModeToUseVector();
250   myGlyph->SetScaleModeToScaleByVector();
251   myGlyph->SetColorModeToColorByScalar();
252 }
253
254
255 //----------------------------------------------------------------------------
256 vtkDataSet* 
257 VISU_VectorsPL
258 ::InsertCustomPL()
259 {
260   return myGlyph->GetOutput();
261 }
262
263
264 //----------------------------------------------------------------------------
265 void
266 VISU_VectorsPL
267 ::Update()
268 {
269   switch (myTypeGlyph) {
270   case ARROW: {
271     myGlyphSource->SetGlyphTypeToArrow();
272     myGlyphSource->SetFilled(0);
273     switch (myPosGlyph) {
274     case TAIL:
275       myGlyphSource->SetCenter(0.5, 0.0, 0.0);
276       break;
277     case HEAD:
278       myGlyphSource->SetCenter(-0.5, 0.0, 0.0);
279       break;
280     case CENTER:
281       myGlyphSource->SetCenter(0.0, 0.0, 0.0);
282     }
283     myGlyph->SetSource(myGlyphSource->GetOutput());
284   }
285     break;
286   case CONE2:
287   case CONE6: {
288     if (myTypeGlyph == CONE2)
289       myConeSource->SetResolution(3);
290     else
291       myConeSource->SetResolution(7);
292     myConeSource->SetHeight(1.0);
293     myConeSource->SetRadius(.1);
294
295     switch (myPosGlyph) {
296     case TAIL:
297       myConeSource->SetCenter(0.5, 0.0, 0.0);
298       break;
299     case HEAD:
300       myConeSource->SetCenter(-0.5, 0.0, 0.0);
301       break;
302     case CENTER:
303       myConeSource->SetCenter(0.0, 0.0, 0.0);
304     }
305     myGlyph->SetSource(myConeSource->GetOutput());
306   }
307     break;
308   case NONE:
309   default: {
310     myGlyph->SetSource(myLineSource->GetOutput());
311   }
312   }
313
314   Superclass::Update();
315 }
316
317
318 //----------------------------------------------------------------------------
319 unsigned long int
320 VISU_VectorsPL
321 ::GetMemorySize()
322 {
323   unsigned long int aSize = Superclass::GetMemorySize();
324
325   if(vtkDataSet* aDataSet = myGlyph->GetOutput())
326     aSize += aDataSet->GetActualMemorySize() * 1024;
327   
328   if(vtkDataSet* aDataSet = myCenters->GetOutput())
329     aSize += aDataSet->GetActualMemorySize() * 1024;
330
331   return aSize;
332 }
333
334
335 //----------------------------------------------------------------------------
336 vtkDataSet* 
337 VISU_VectorsPL
338 ::GetOutput()
339 {
340   myGlyph->Update();
341
342   return myGlyph->GetOutput();
343 }
344
345
346 //----------------------------------------------------------------------------
347 void
348 VISU_VectorsPL
349 ::SetMapScale(vtkFloatingPointType theMapScale)
350 {
351   VISU_ScalarMapPL::SetMapScale(theMapScale);
352   myMapScaleFactor = theMapScale;
353
354   myGlyph->SetScaleFactor( myScaleFactor*theMapScale );
355
356   Modified();
357 }
358
359
360 //----------------------------------------------------------------------------