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