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