Salome HOME
7696d62e5d5453ae344ec30efd6e543808ef9011
[modules/visu.git] / src / PIPELINE / VISU_Plot3DPL.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_Plot3DPL.hxx"
29 #include "VISU_CutPlanesPL.hxx"
30 #include "VISU_PipeLineUtils.hxx"
31
32 #include <vtkAppendPolyData.h>
33 #include <vtkCutter.h>
34 #include <vtkPlane.h>
35
36 #include <vtkCellDataToPointData.h>
37 #include <vtkGeometryFilter.h>
38 #include <vtkContourFilter.h>
39 #include <vtkWarpScalar.h>
40 #include <vtkOutlineFilter.h>
41
42
43 //----------------------------------------------------------------------------
44 vtkStandardNewMacro(VISU_Plot3DPL);
45
46
47 //----------------------------------------------------------------------------
48 VISU_Plot3DPL
49 ::VISU_Plot3DPL():
50   myCellDataToPointData(vtkCellDataToPointData::New()),
51   myAppendPolyData(vtkAppendPolyData::New()),
52   myGeometryFilter(vtkGeometryFilter::New()),
53   myContourFilter(vtkContourFilter::New()),
54   myWarpScalar(vtkWarpScalar::New()),
55   myOrientation(VISU_CutPlanesPL::YZ),
56   myIsRelative(true),
57   myIsContour(false),
58   myPosition(0.5),
59   myScaleFactor(1.0)
60 {
61   SetIsShrinkable(false);
62
63   myCellDataToPointData->Delete();
64   myAppendPolyData->Delete();
65   myGeometryFilter->Delete();
66   myContourFilter->Delete();
67   myWarpScalar->Delete();
68
69   myAngle[0] = myAngle[1] = myAngle[2] = 0.0;
70
71   SetNumberOfContours(32);
72 }
73
74
75 //----------------------------------------------------------------------------
76 VISU_Plot3DPL
77 ::~VISU_Plot3DPL()
78 {}
79
80
81 //----------------------------------------------------------------------------
82 unsigned long int 
83 VISU_Plot3DPL
84 ::GetMTime()
85 {
86   unsigned long int aTime = Superclass::GetMTime();
87
88   aTime = std::max(aTime, myCellDataToPointData->GetMTime());
89   aTime = std::max(aTime, myAppendPolyData->GetMTime());
90   aTime = std::max(aTime, myGeometryFilter->GetMTime());
91   aTime = std::max(aTime, myContourFilter->GetMTime());
92   aTime = std::max(aTime, myWarpScalar->GetMTime());
93
94   return aTime;
95 }
96
97
98 //----------------------------------------------------------------------------
99 void
100 VISU_Plot3DPL
101 ::DoShallowCopy(VISU_PipeLine *thePipeLine,
102                 bool theIsCopyInput)
103 {
104   Superclass::DoShallowCopy(thePipeLine, theIsCopyInput);
105
106   if(VISU_Plot3DPL *aPipeLine = dynamic_cast<VISU_Plot3DPL*>(thePipeLine)){
107     SetOrientation (aPipeLine->GetPlaneOrientation(),
108                     aPipeLine->GetRotateX(), aPipeLine->GetRotateY());
109     SetPlanePosition (aPipeLine->GetPlanePosition(),
110                       aPipeLine->IsPositionRelative() );
111     SetScaleFactor( aPipeLine->GetScaleFactor() );
112     SetContourPrs( aPipeLine->GetIsContourPrs() );
113     SetNumberOfContours( aPipeLine->GetNumberOfContours() );
114   }
115 }
116
117
118 //----------------------------------------------------------------------------
119 VISU_CutPlanesPL::PlaneOrientation
120 VISU_Plot3DPL
121 ::GetOrientation(vtkDataSet* theDataSet)
122 {
123   theDataSet->Update();
124
125   vtkFloatingPointType aBounds[6];
126   theDataSet->GetBounds(aBounds);
127   vtkFloatingPointType aDelta[3] = {aBounds[1] - aBounds[0], aBounds[3] - aBounds[2], aBounds[5] - aBounds[4]};
128
129   if(aDelta[0] >= aDelta[1] && aDelta[0] >= aDelta[2])
130     if(aDelta[1] >= aDelta[2])
131       return VISU_CutPlanesPL::XY;
132     else
133       return VISU_CutPlanesPL::ZX;
134
135   if(aDelta[1] >= aDelta[0] && aDelta[1] >= aDelta[2])
136     if(aDelta[0] >= aDelta[2])
137       return VISU_CutPlanesPL::XY;
138     else
139       return VISU_CutPlanesPL::YZ;
140
141   if(aDelta[2] >= aDelta[0] && aDelta[2] >= aDelta[1])
142     if(aDelta[0] >= aDelta[1])
143       return VISU_CutPlanesPL::ZX;
144     else
145       return VISU_CutPlanesPL::YZ;
146
147   return VISU_CutPlanesPL::XY;
148 }
149
150
151 //----------------------------------------------------------------------------
152 vtkFloatingPointType
153 VISU_Plot3DPL
154 ::GetScaleFactor(vtkDataSet* theDataSet)
155 {
156   theDataSet->Update();
157   vtkFloatingPointType aLength = theDataSet->GetLength(); // diagonal length
158
159   vtkFloatingPointType aScalarRange[2];
160   GetSourceRange(aScalarRange);
161
162   static vtkFloatingPointType EPS = 0.3;
163   vtkFloatingPointType aRange = aScalarRange[1];
164   if(aRange > 0.0)
165     return aLength / aRange * EPS;
166
167   return 0.0;
168 }
169
170
171 //----------------------------------------------------------------------------
172 void
173 VISU_Plot3DPL
174 ::Init()
175 {
176   Superclass::Init();
177
178   myOrientation = GetOrientation(GetMergedInput());
179   SetScaleFactor(GetScaleFactor(GetMergedInput()));
180 }
181
182
183 //----------------------------------------------------------------------------
184 vtkDataSet*
185 VISU_Plot3DPL
186 ::InsertCustomPL()
187 {
188   return myAppendPolyData->GetOutput();
189 }
190
191
192 //----------------------------------------------------------------------------
193 void
194 VISU_Plot3DPL
195 ::Update()
196 {
197   vtkFloatingPointType aPlaneNormal[3];
198   vtkFloatingPointType anOrigin[3];
199   GetBasePlane( anOrigin, aPlaneNormal );
200
201   vtkPolyData* aPolyData = 0;
202   vtkCutter *aCutPlane = 0;
203   vtkDataSet* aDataSet = GetMergedInput();
204
205   if ( !IsPlanarInput() )
206   {
207     aCutPlane = vtkCutter::New();
208     aCutPlane->SetInput(aDataSet);
209
210     vtkPlane *aPlane = vtkPlane::New();
211     aPlane->SetOrigin(anOrigin);
212     aPlane->SetNormal(aPlaneNormal);
213
214     aCutPlane->SetCutFunction(aPlane);
215     aPlane->Delete();
216
217     aPolyData = aCutPlane->GetOutput();
218     aPolyData->Update();
219   }
220
221   if ( !aPolyData || aPolyData->GetNumberOfCells() == 0 ) {
222     myGeometryFilter->SetInput(aDataSet);
223     aPolyData = myGeometryFilter->GetOutput();
224     aPolyData->Update();
225   }
226   if ( !myIsContour ) // surface prs
227   {
228     if(VISU::IsDataOnCells(aPolyData)) {
229       myCellDataToPointData->SetInput(aPolyData);
230       myCellDataToPointData->PassCellDataOn();
231       myWarpScalar->SetInput(myCellDataToPointData->GetPolyDataOutput());
232     }else
233       myWarpScalar->SetInput(aPolyData);
234   }
235   else // contour prs
236   {
237     if(VISU::IsDataOnCells(aPolyData)) {
238       myCellDataToPointData->SetInput(aPolyData);
239       myCellDataToPointData->PassCellDataOn();
240       myContourFilter->SetInput(myCellDataToPointData->GetOutput());
241     }else
242       myContourFilter->SetInput(aPolyData);
243
244     vtkFloatingPointType aScalarRange[2];
245     GetSourceRange(aScalarRange);
246
247     myContourFilter->GenerateValues(GetNumberOfContours(),aScalarRange);
248     myWarpScalar->SetInput(myContourFilter->GetOutput());
249   }
250
251   VISU_CutPlanesPL::ClearAppendPolyData(myAppendPolyData.GetPointer());
252   myAppendPolyData->AddInput(myWarpScalar->GetPolyDataOutput());
253
254   if ( aCutPlane )
255     aCutPlane->Delete();
256
257   myWarpScalar->SetNormal(aPlaneNormal);
258
259   Superclass::Update();
260 }
261
262
263 //----------------------------------------------------------------------------
264 unsigned long int
265 VISU_Plot3DPL
266 ::GetMemorySize()
267 {
268   unsigned long int aSize = Superclass::GetMemorySize();
269
270   if(vtkDataObject* aDataObject = myGeometryFilter->GetInput())
271     aSize += aDataObject->GetActualMemorySize() * 1024;
272   
273   if(myCellDataToPointData->GetInput())
274     if(vtkDataSet* aDataSet = myCellDataToPointData->GetOutput())
275       aSize += aDataSet->GetActualMemorySize() * 1024;
276
277   if(vtkDataObject* aDataObject = myContourFilter->GetInput())
278     aSize += aDataObject->GetActualMemorySize() * 1024;
279
280   if(vtkDataObject* aDataObject = myWarpScalar->GetInput())
281     aSize += aDataObject->GetActualMemorySize() * 1024;
282
283   int anEnd = myAppendPolyData->GetNumberOfInputConnections(0);
284   for(int anId = 0; anId < anEnd; anId++){
285     if(vtkDataObject* aDataObject = myAppendPolyData->GetInput(anId))
286       aSize += aDataObject->GetActualMemorySize() * 1024;
287   }
288
289   return aSize;
290 }
291
292
293 //----------------------------------------------------------------------------
294 void
295 VISU_Plot3DPL
296 ::SetNumberOfContours(int theNumber)
297 {
298   myContourFilter->SetNumberOfContours(theNumber);
299 }
300
301
302 //----------------------------------------------------------------------------
303 int
304 VISU_Plot3DPL
305 ::GetNumberOfContours()
306 {
307   return myContourFilter->GetNumberOfContours();
308 }
309
310
311 //----------------------------------------------------------------------------
312 void
313 VISU_Plot3DPL
314 ::SetScaleFactor(vtkFloatingPointType theScaleFactor)
315 {
316   myScaleFactor = theScaleFactor;
317   myWarpScalar->SetScaleFactor(theScaleFactor);
318 }
319
320
321 //----------------------------------------------------------------------------
322 vtkFloatingPointType
323 VISU_Plot3DPL
324 ::GetScaleFactor()
325 {
326   return myScaleFactor;
327 }
328
329
330 //----------------------------------------------------------------------------
331 void
332 VISU_Plot3DPL::
333 SetContourPrs(bool theIsContourPrs )
334 {
335   if(myIsContour == theIsContourPrs)
336     return;
337
338   myIsContour = theIsContourPrs;
339   Modified();
340 }
341
342
343 //----------------------------------------------------------------------------
344 bool
345 VISU_Plot3DPL
346 ::GetIsContourPrs()
347 {
348   return myIsContour;
349 }
350
351
352 //----------------------------------------------------------------------------
353 void
354 VISU_Plot3DPL
355 ::SetPlanePosition(vtkFloatingPointType thePosition,
356                    bool theIsRelative)
357 {
358   bool anIsSameValue = VISU::CheckIsSameValue(myIsRelative, theIsRelative);
359   anIsSameValue &= (myPosition == thePosition);
360   if(anIsSameValue)
361     return;
362
363   myIsRelative = theIsRelative;
364   myPosition = thePosition;
365   Modified();
366 }
367
368
369 //----------------------------------------------------------------------------
370 bool
371 VISU_Plot3DPL
372 ::IsPositionRelative()
373 {
374   return myIsRelative;
375 }
376
377
378 //----------------------------------------------------------------------------
379 VISU_CutPlanesPL::PlaneOrientation
380 VISU_Plot3DPL
381 ::GetPlaneOrientation()
382 {
383   return myOrientation;
384 }
385
386
387 //----------------------------------------------------------------------------
388 vtkFloatingPointType
389 VISU_Plot3DPL::
390 GetRotateX()
391 {
392   switch(myOrientation){
393   case VISU_CutPlanesPL::XY: return myAngle[0];
394   case VISU_CutPlanesPL::YZ: return myAngle[1];
395   case VISU_CutPlanesPL::ZX: return myAngle[2];
396   }
397   return 0;
398 }
399
400
401 //----------------------------------------------------------------------------
402 vtkFloatingPointType
403 VISU_Plot3DPL::
404 GetRotateY(){
405   switch(myOrientation){
406   case VISU_CutPlanesPL::XY: return myAngle[1];
407   case VISU_CutPlanesPL::YZ: return myAngle[2];
408   case VISU_CutPlanesPL::ZX: return myAngle[0];
409   }
410   return 0;
411 }
412
413
414 //----------------------------------------------------------------------------
415 void
416 VISU_Plot3DPL::
417 SetOrientation(VISU_CutPlanesPL::PlaneOrientation theOrientation,
418                vtkFloatingPointType theXAngle,
419                vtkFloatingPointType theYAngle)
420 {
421   bool anIsSameValue = VISU::CheckIsSameValue(GetRotateX(), theXAngle);
422   anIsSameValue &= VISU::CheckIsSameValue(GetRotateY(), theYAngle);
423   anIsSameValue &= (myOrientation == theOrientation);
424   if(anIsSameValue)
425     return;
426
427   switch(theOrientation){
428   case VISU_CutPlanesPL::XY: myAngle[0] = theXAngle; break;
429   case VISU_CutPlanesPL::YZ: myAngle[1] = theXAngle; break;
430   case VISU_CutPlanesPL::ZX: myAngle[2] = theXAngle; break;
431   }
432
433   switch(theOrientation){
434   case VISU_CutPlanesPL::XY: myAngle[1] = theYAngle; break;
435   case VISU_CutPlanesPL::YZ: myAngle[2] = theYAngle; break;
436   case VISU_CutPlanesPL::ZX: myAngle[0] = theYAngle; break;
437   }
438
439   myOrientation = theOrientation;
440   Modified();
441 }
442
443
444 //----------------------------------------------------------------------------
445 vtkFloatingPointType
446 VISU_Plot3DPL
447 ::GetPlanePosition()
448 {
449   return myPosition;
450 }
451
452 //=======================================================================
453 //function : GetBasePlane
454 //purpose  :
455 //=======================================================================
456 void
457 VISU_Plot3DPL
458 ::GetBasePlane(vtkFloatingPointType theOrigin[3],
459                vtkFloatingPointType theNormal[3],
460                bool  theCenterOrigine )
461 {
462   VISU_CutPlanesPL::GetDir(theNormal,myAngle,myOrientation);
463
464   vtkFloatingPointType aPosition = myPosition;
465   vtkFloatingPointType aBounds[6], aBoundPrj[3];
466   if ( myIsRelative )
467   {
468     GetInput()->GetBounds(aBounds);
469     VISU_CutPlanesPL::GetBoundProject(aBoundPrj,aBounds,theNormal);
470     aPosition = aBoundPrj[0] + aBoundPrj[2]*myPosition;
471   }
472   VISU::Mul(theNormal,aPosition,theOrigin);
473
474   if ( theCenterOrigine ) {
475     // move theOrigin to the center of aBounds projections to the plane
476     GetMergedInput()->GetBounds(aBounds);
477     vtkFloatingPointType boundPoints[8][3] = {
478       {aBounds[0],aBounds[2],aBounds[4]},
479       {aBounds[1],aBounds[2],aBounds[4]},
480       {aBounds[0],aBounds[3],aBounds[4]},
481       {aBounds[1],aBounds[3],aBounds[4]},
482       {aBounds[0],aBounds[2],aBounds[5]},
483       {aBounds[1],aBounds[2],aBounds[5]},
484       {aBounds[0],aBounds[3],aBounds[5]},
485       {aBounds[1],aBounds[3],aBounds[5]}};
486     vtkFloatingPointType newOrigin[3] = { 0,0,0 };
487     for(int i = 0; i < 8; i++) {
488       vtkFloatingPointType proj[3];
489       vtkPlane::ProjectPoint( boundPoints[i], theOrigin, theNormal, proj );
490       newOrigin[0] += proj[0];
491       newOrigin[1] += proj[1];
492       newOrigin[2] += proj[2];
493     }
494     theOrigin[0] = newOrigin[0] / 8.;
495     theOrigin[1] = newOrigin[1] / 8.;
496     theOrigin[2] = newOrigin[2] / 8.;
497   }
498 }
499
500 //=======================================================================
501 //function : GetMinMaxPosition
502 //purpose  : return absolute position range
503 //=======================================================================
504 void
505 VISU_Plot3DPL
506 ::GetMinMaxPosition( vtkFloatingPointType& minPos, 
507                      vtkFloatingPointType& maxPos )
508 {
509   vtkFloatingPointType aBounds[6], aBoundPrj[3], aNormal[3];
510   VISU_CutPlanesPL::GetDir(aNormal,myAngle,myOrientation);
511   GetInput()->GetBounds(aBounds);
512   VISU_CutPlanesPL::GetBoundProject(aBoundPrj,aBounds,aNormal);
513   minPos = aBoundPrj[0];
514   maxPos = aBoundPrj[1];
515 }
516
517 //=======================================================================
518 //function : SetMapScale
519 //purpose  :
520 //=======================================================================
521
522 void 
523 VISU_Plot3DPL
524 ::SetMapScale(vtkFloatingPointType theMapScale)
525 {
526   Superclass::SetMapScale(theMapScale);
527
528   if ( myIsContour ) {
529     vtkFloatingPointType aRange[2];
530     GetSourceRange(aRange);
531     vtkFloatingPointType aNewRange[] = { aRange[1] - theMapScale*(aRange[1]-aRange[0]), aRange[1] };
532     myContourFilter->GenerateValues(GetNumberOfContours(),aNewRange);
533   }
534   myWarpScalar->SetScaleFactor(myScaleFactor*theMapScale);
535
536   Modified();
537 }