Salome HOME
Fix for Bug PAL8597:
[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.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org
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 using namespace std;
43
44 vtkStandardNewMacro(VISU_Plot3DPL);
45
46 VISU_Plot3DPL::VISU_Plot3DPL():
47   myCellDataToPointData(vtkCellDataToPointData::New(),true),
48   myAppendPolyData(vtkAppendPolyData::New(),true),
49   myGeometryFilter(vtkGeometryFilter::New(),true),
50   myContourFilter(vtkContourFilter::New(),true),
51   myWarpScalar(vtkWarpScalar::New(),true),
52   myOrientation(VISU_CutPlanesPL::YZ),
53   myIsRelative(true),
54   myIsContour(false),
55   myPosition(0.5),
56   myScaleFactor(1.)
57 {
58   myAngle[0] = myAngle[1] = myAngle[2] = 0.;
59   SetNumberOfContours(32);
60   myIsShrinkable = false;
61 }
62
63 VISU_Plot3DPL::~VISU_Plot3DPL()
64 {
65 }
66
67 void
68 VISU_Plot3DPL::
69 ShallowCopy(VISU_PipeLine *thePipeLine)
70 {
71   VISU_ScalarMapPL::ShallowCopy(thePipeLine);
72   if(VISU_Plot3DPL *aPipeLine = dynamic_cast<VISU_Plot3DPL*>(thePipeLine)){
73     SetOrientation (aPipeLine->GetPlaneOrientation(),
74                     aPipeLine->GetRotateX(), aPipeLine->GetRotateY());
75     SetPlanePosition (aPipeLine->GetPlanePosition(),
76                       aPipeLine->IsPositionRelative() );
77     SetScaleFactor( aPipeLine->GetScaleFactor() );
78     SetContourPrs( aPipeLine->GetIsContourPrs() );
79     SetNumberOfContours( aPipeLine->GetNumberOfContours() );
80   }
81 }
82
83 VISU_CutPlanesPL::PlaneOrientation
84 VISU_Plot3DPL::
85 GetOrientation(vtkDataSet* theDataSet)
86 {
87   theDataSet->Update();
88
89   float aBounds[6];
90   theDataSet->GetBounds(aBounds);
91   float aDelta[3] = {aBounds[1] - aBounds[0], aBounds[3] - aBounds[2], aBounds[5] - aBounds[4]};
92
93   if(aDelta[0] >= aDelta[1] && aDelta[0] >= aDelta[2])
94     if(aDelta[1] >= aDelta[2])
95       return VISU_CutPlanesPL::XY;
96     else
97       return VISU_CutPlanesPL::ZX;
98
99   if(aDelta[1] >= aDelta[0] && aDelta[1] >= aDelta[2])
100     if(aDelta[0] >= aDelta[2])
101       return VISU_CutPlanesPL::XY;
102     else
103       return VISU_CutPlanesPL::YZ;
104
105   if(aDelta[2] >= aDelta[0] && aDelta[2] >= aDelta[1])
106     if(aDelta[0] >= aDelta[1])
107       return VISU_CutPlanesPL::ZX;
108     else
109       return VISU_CutPlanesPL::YZ;
110
111   return VISU_CutPlanesPL::XY;
112 }
113
114 float
115 VISU_Plot3DPL::
116 GetScaleFactor(vtkDataSet* theDataSet)
117 {
118   theDataSet->Update();
119   float aLength = theDataSet->GetLength(); // diagonal length
120
121   float aScalarRange[2];
122   GetSourceRange(aScalarRange);
123
124   static float EPS = 0.3;
125   float aRange = aScalarRange[1];
126   if(aRange > 0.0)
127     return aLength / aRange * EPS;
128
129   return 0.0;
130 }
131
132 void
133 VISU_Plot3DPL::
134 Init()
135 {
136   VISU_ScalarMapPL::Init();
137
138   myOrientation = GetOrientation(GetInput2());
139   SetScaleFactor(GetScaleFactor(GetInput2()));
140 }
141
142 VISU_ScalarMapPL::THook*
143 VISU_Plot3DPL::
144 DoHook()
145 {
146   return myAppendPolyData->GetOutput();
147 }
148
149 void
150 VISU_Plot3DPL::
151 Update()
152 {
153   float aPlaneNormal[3];
154   float anOrigin[3];
155   GetBasePlane( anOrigin, aPlaneNormal );
156
157   vtkPolyData* aPolyData = 0;
158   vtkCutter *aCutPlane = 0;
159   vtkUnstructuredGrid* anUnstructuredGrid =
160     myFieldTransform->GetUnstructuredGridOutput();
161
162   if ( !IsPlanarInput() )
163   {
164     aCutPlane = vtkCutter::New();
165     aCutPlane->SetInput(anUnstructuredGrid);
166
167     vtkPlane *aPlane = vtkPlane::New();
168     aPlane->SetOrigin(anOrigin);
169     aPlane->SetNormal(aPlaneNormal);
170
171     aCutPlane->SetCutFunction(aPlane);
172     aPlane->Delete();
173
174     aPolyData = aCutPlane->GetOutput();
175     aPolyData->Update();
176   }
177
178   if ( !aPolyData || aPolyData->GetNumberOfCells() == 0 ) {
179     myGeometryFilter->SetInput(anUnstructuredGrid);
180     aPolyData = myGeometryFilter->GetOutput();
181     aPolyData->Update();
182   }
183   if ( !myIsContour ) // surface prs
184   {
185     if(aPolyData->GetCellData()->GetNumberOfArrays()) {
186       myCellDataToPointData->SetInput(aPolyData);
187       myCellDataToPointData->PassCellDataOn();
188       myWarpScalar->SetInput(myCellDataToPointData->GetPolyDataOutput());
189     }else
190       myWarpScalar->SetInput(aPolyData);
191   }
192   else // contour prs
193   {
194     if(aPolyData->GetCellData()->GetNumberOfArrays()) {
195       myCellDataToPointData->SetInput(aPolyData);
196       myCellDataToPointData->PassCellDataOn();
197       myContourFilter->SetInput(myCellDataToPointData->GetOutput());
198     }else
199       myContourFilter->SetInput(aPolyData);
200
201     float aScalarRange[2];
202     GetSourceRange(aScalarRange);
203
204     myContourFilter->GenerateValues(GetNumberOfContours(),aScalarRange);
205     myWarpScalar->SetInput(myContourFilter->GetOutput());
206   }
207
208   VISU_CutPlanesPL::ClearAppendPolyData(myAppendPolyData);
209   myAppendPolyData->AddInput(myWarpScalar->GetPolyDataOutput());
210
211   if ( aCutPlane )
212     aCutPlane->Delete();
213
214   myWarpScalar->SetNormal(aPlaneNormal);
215
216   VISU_ScalarMapPL::Update();
217 }
218
219 void
220 VISU_Plot3DPL::
221 SetNumberOfContours(int theNumber)
222 {
223   myContourFilter->SetNumberOfContours(theNumber);
224 }
225
226 int
227 VISU_Plot3DPL::
228 GetNumberOfContours() const
229 {
230   return myContourFilter->GetNumberOfContours();
231 }
232
233 void
234 VISU_Plot3DPL::
235 SetScaleFactor(float theScaleFactor)
236 {
237   myScaleFactor = theScaleFactor;
238   myWarpScalar->SetScaleFactor(theScaleFactor);
239 }
240
241 float
242 VISU_Plot3DPL::
243 GetScaleFactor() const
244 {
245   return myScaleFactor;
246 }
247
248 void
249 VISU_Plot3DPL::
250 SetPlanePosition(float thePosition,
251                  bool theIsRelative)
252 {
253   myIsRelative = theIsRelative;
254   myPosition = thePosition;
255 }
256
257 bool
258 VISU_Plot3DPL::
259 IsPositionRelative() const
260 {
261   return myIsRelative;
262 }
263
264 VISU_CutPlanesPL::PlaneOrientation
265 VISU_Plot3DPL::
266 GetPlaneOrientation() const
267 {
268   return myOrientation;
269 }
270
271
272 float
273 VISU_Plot3DPL::
274 GetRotateX()
275 {
276   switch(myOrientation){
277   case VISU_CutPlanesPL::XY: return myAngle[0];
278   case VISU_CutPlanesPL::YZ: return myAngle[1];
279   case VISU_CutPlanesPL::ZX: return myAngle[2];
280   }
281   return 0;
282 }
283
284 float
285 VISU_Plot3DPL::
286 GetRotateY(){
287   switch(myOrientation){
288   case VISU_CutPlanesPL::XY: return myAngle[1];
289   case VISU_CutPlanesPL::YZ: return myAngle[2];
290   case VISU_CutPlanesPL::ZX: return myAngle[0];
291   }
292   return 0;
293 }
294
295 void
296 VISU_Plot3DPL::
297 SetOrientation(VISU_CutPlanesPL::PlaneOrientation theOrientation,
298                float theXAngle,
299                float theYAngle)
300 {
301   switch(theOrientation){
302   case VISU_CutPlanesPL::XY: myAngle[0] = theXAngle; break;
303   case VISU_CutPlanesPL::YZ: myAngle[1] = theXAngle; break;
304   case VISU_CutPlanesPL::ZX: myAngle[2] = theXAngle; break;
305   }
306   switch(theOrientation){
307   case VISU_CutPlanesPL::XY: myAngle[1] = theYAngle; break;
308   case VISU_CutPlanesPL::YZ: myAngle[2] = theYAngle; break;
309   case VISU_CutPlanesPL::ZX: myAngle[0] = theYAngle; break;
310   }
311   myOrientation = theOrientation;
312 }
313
314 float
315 VISU_Plot3DPL::
316 GetPlanePosition() const
317 {
318   return myPosition;
319 }
320
321 //=======================================================================
322 //function : GetBasePlane
323 //purpose  :
324 //=======================================================================
325
326 void VISU_Plot3DPL::GetBasePlane(float theOrigin[3],
327                                  float theNormal[3],
328                                  bool  theCenterOrigine ) const
329 {
330   VISU_CutPlanesPL::GetDir(theNormal,myAngle,myOrientation);
331
332   float aPosition = myPosition;
333   float aBounds[6], aBoundPrj[3];
334   if ( myIsRelative )
335   {
336     GetInput()->GetBounds(aBounds);
337     VISU_CutPlanesPL::GetBoundProject(aBoundPrj,aBounds,theNormal);
338     aPosition = aBoundPrj[0] + aBoundPrj[2]*myPosition;
339   }
340   VISU::Mul(theNormal,aPosition,theOrigin);
341
342   if ( theCenterOrigine ) {
343     // move theOrigin to the center of aBounds projections to the plane
344     GetInput2()->GetBounds(aBounds);
345     float boundPoints[8][3] = {
346       {aBounds[0],aBounds[2],aBounds[4]},
347       {aBounds[1],aBounds[2],aBounds[4]},
348       {aBounds[0],aBounds[3],aBounds[4]},
349       {aBounds[1],aBounds[3],aBounds[4]},
350       {aBounds[0],aBounds[2],aBounds[5]},
351       {aBounds[1],aBounds[2],aBounds[5]},
352       {aBounds[0],aBounds[3],aBounds[5]},
353       {aBounds[1],aBounds[3],aBounds[5]}};
354     float newOrigin[3] = { 0,0,0 };
355     for(int i = 0; i < 8; i++) {
356       float proj[3];
357       vtkPlane::ProjectPoint( boundPoints[i], theOrigin, theNormal, proj );
358       newOrigin[0] += proj[0];
359       newOrigin[1] += proj[1];
360       newOrigin[2] += proj[2];
361     }
362     theOrigin[0] = newOrigin[0] / 8.;
363     theOrigin[1] = newOrigin[1] / 8.;
364     theOrigin[2] = newOrigin[2] / 8.;
365   }
366 }
367
368 //=======================================================================
369 //function : GetMinMaxPosition
370 //purpose  : return absolute position range
371 //=======================================================================
372
373 void VISU_Plot3DPL::GetMinMaxPosition( float& minPos, float& maxPos ) const
374 {
375   float aBounds[6], aBoundPrj[3], aNormal[3];
376   VISU_CutPlanesPL::GetDir(aNormal,myAngle,myOrientation);
377   GetInput()->GetBounds(aBounds);
378   VISU_CutPlanesPL::GetBoundProject(aBoundPrj,aBounds,aNormal);
379   minPos = aBoundPrj[0];
380   maxPos = aBoundPrj[1];
381 }
382
383 //=======================================================================
384 //function : SetMapScale
385 //purpose  :
386 //=======================================================================
387
388 void VISU_Plot3DPL::SetMapScale(float theMapScale)
389 {
390   VISU_ScalarMapPL::SetMapScale(theMapScale);
391
392   if ( myIsContour ) {
393     float aRange[2];
394     GetSourceRange(aRange);
395     float aNewRange[] = { aRange[1] - theMapScale*(aRange[1]-aRange[0]), aRange[1] };
396     myContourFilter->GenerateValues(GetNumberOfContours(),aNewRange);
397   }
398   myWarpScalar->SetScaleFactor(myScaleFactor*theMapScale);
399
400   Modified();
401 }