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