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