Salome HOME
63dd3e5a42966a4dee8b3c794f1fe08e2b8a42a2
[modules/visu.git] / src / PIPELINE / VISU_CutPlanesPL.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_CutPlanesPL.hxx"
29 #include "VISU_PipeLineUtils.hxx"
30 #include "VTKViewer_GeometryFilter.h"
31
32 #include <vtkAppendPolyData.h>
33 #include <vtkCutter.h>
34 #include <vtkPlane.h>
35
36 using namespace std;
37
38 static float EPS = 1.0E-3;
39
40 vtkStandardNewMacro(VISU_CutPlanesPL);
41
42 VISU_CutPlanesPL::VISU_CutPlanesPL(){
43   myAppendPolyData = vtkAppendPolyData::New();
44   myIsShrinkable = false;
45 }
46
47 VISU_CutPlanesPL::~VISU_CutPlanesPL(){
48   myAppendPolyData->Delete();
49 }
50
51 void VISU_CutPlanesPL::ShallowCopy(VISU_PipeLine *thePipeLine){
52   if(VISU_CutPlanesPL *aPipeLine = dynamic_cast<VISU_CutPlanesPL*>(thePipeLine)){
53     SetOrientation(aPipeLine->GetPlaneOrientation(),
54                    aPipeLine->GetRotateX(),aPipeLine->GetRotateY());
55     SetDisplacement(aPipeLine->GetDisplacement());
56     SetNbParts(aPipeLine->GetNbParts());
57     for (int i = 0, iend = GetNbParts(); i < iend; i++)
58       if(!aPipeLine->IsPartDefault(i))  SetPartPosition(i, aPipeLine->GetPartPosition(i));
59   }
60   VISU_ScalarMapPL::ShallowCopy(thePipeLine);
61 }
62
63 void VISU_CutPlanesPL::Init(){
64   VISU_ScalarMapPL::Init();
65
66   SetNbParts(10);
67   myBasePlane[0] = YZ;
68   myDisplacement[0] = 0.5;
69   myAng[0][0] = myAng[0][1] = myAng[0][2] = 0.0;
70 }
71
72 VISU_ScalarMapPL::THook* VISU_CutPlanesPL::DoHook(){
73   return myAppendPolyData->GetOutput();
74 }
75
76 void VISU_CutPlanesPL::Update(){
77   ClearAppendPolyData(myAppendPolyData);
78   SetPartPosition();
79   float aDir[3];
80   GetDir(aDir,myAng[0],myBasePlane[0]);
81   float aBounds[6];
82   GetInput2()->GetBounds(aBounds);
83   vtkDataSet* aDataSet = myFieldTransform->GetUnstructuredGridOutput();
84   CutWithPlanes(myAppendPolyData,aDataSet,myNbParts,aDir,aBounds,
85                 myPartPosition,myPartCondition,myDisplacement[0]);
86
87   VISU_ScalarMapPL::Update();
88 }
89
90 void VISU_CutPlanesPL::SetPartPosition(int theNum){
91   for(int i = 0; i < myNbParts; i++)
92     myPartPosition[i] = GetPartPosition(i,theNum);
93 }
94
95 void VISU_CutPlanesPL::ClearAppendPolyData(vtkAppendPolyData *theAppendPolyData){
96   int iEnd = theAppendPolyData->GetNumberOfInputs();
97   for(int i = iEnd-1; i >= 0; i--)
98     theAppendPolyData->RemoveInput(theAppendPolyData->GetInput(i));
99 }
100
101 float* VISU_CutPlanesPL::GetRx(float theRx[3][3], float thaAng){
102   theRx[0][0] = 1.0;            theRx[0][1] = 0.0;            theRx[0][2] = 0.0;
103   theRx[1][0] = 0.0;            theRx[1][1] = cos(thaAng);    theRx[1][2] = -sin(thaAng);
104   theRx[2][0] = 0.0;            theRx[2][1] = sin(thaAng);    theRx[2][2] = cos(thaAng);
105   return theRx[0];
106 }
107
108
109 float* VISU_CutPlanesPL::GetRy(float theRy[3][3], float thaAng){
110   theRy[0][0] = cos(thaAng);    theRy[0][1] = 0.0;            theRy[0][2] = sin(thaAng);
111   theRy[1][0] = 0.0;            theRy[1][1] = 1.0;            theRy[1][2] = 0.0;
112   theRy[2][0] = -sin(thaAng);   theRy[2][1] = 0.0;            theRy[2][2] = cos(thaAng);
113   return theRy[0];
114 }
115
116
117 float* VISU_CutPlanesPL::GetRz(float theRz[3][3], float thaAng){
118   theRz[0][0] = cos(thaAng);    theRz[0][1] = -sin(thaAng);   theRz[0][2] = 0.0;
119   theRz[1][0] = sin(thaAng);    theRz[1][1] = cos(thaAng);    theRz[1][2] = 0.0;
120   theRz[2][0] = 0.0;            theRz[2][1] = 0.0;            theRz[2][2] = 1.0;
121   return theRz[0];
122 }
123
124
125 void VISU_CutPlanesPL::CorrectPnt(float thePnt[3], const float BoundPrj[6]){
126   for(int i = 0, j = 0; i < 3; ++i, j=2*i){
127     if(thePnt[i] < BoundPrj[j]) thePnt[i] = BoundPrj[j];
128     if(thePnt[i] > BoundPrj[j+1]) thePnt[i] = BoundPrj[j+1];
129   }
130 }
131
132 void VISU_CutPlanesPL::GetBoundProject(float BoundPrj[3], const float BoundBox[6], const float Dir[3]){
133   float BoundPoints[8][3] = { {BoundBox[0],BoundBox[2],BoundBox[4]},
134                               {BoundBox[1],BoundBox[2],BoundBox[4]},
135                               {BoundBox[0],BoundBox[3],BoundBox[4]},
136                               {BoundBox[1],BoundBox[3],BoundBox[4]},
137                               {BoundBox[0],BoundBox[2],BoundBox[5]},
138                               {BoundBox[1],BoundBox[2],BoundBox[5]},
139                               {BoundBox[0],BoundBox[3],BoundBox[5]},
140                               {BoundBox[1],BoundBox[3],BoundBox[5]}};
141   BoundPrj[0] = vtkMath::Dot(Dir,BoundPoints[0]), BoundPrj[1] = BoundPrj[0];
142   for(int i = 1; i < 8; i++){
143     float tmp = vtkMath::Dot(Dir,BoundPoints[i]);
144     if(BoundPrj[1] < tmp) BoundPrj[1] = tmp;
145     if(BoundPrj[0] > tmp) BoundPrj[0] = tmp;
146   }
147   BoundPrj[2] = BoundPrj[1] - BoundPrj[0];
148   BoundPrj[1] = BoundPrj[0] + (1.0 - EPS)*BoundPrj[2];
149   BoundPrj[0] = BoundPrj[0] + EPS*BoundPrj[2];
150   BoundPrj[2] = BoundPrj[1] - BoundPrj[0];
151 }
152
153
154 void VISU_CutPlanesPL::SetOrientation(const VISU_CutPlanesPL::PlaneOrientation& theOrient,
155                                       float theXAng, float theYAng, int theNum)
156 {
157   myBasePlane[theNum] = theOrient;
158   switch(myBasePlane[theNum]){
159   case XY: myAng[theNum][0] = theXAng; break;
160   case YZ: myAng[theNum][1] = theXAng; break;
161   case ZX: myAng[theNum][2] = theXAng; break;
162   }
163   switch(myBasePlane[theNum]){
164   case XY: myAng[theNum][1] = theYAng; break;
165   case YZ: myAng[theNum][2] = theYAng; break;
166   case ZX: myAng[theNum][0] = theYAng; break;
167   }
168 }
169
170
171 const VISU_CutPlanesPL::PlaneOrientation& VISU_CutPlanesPL::GetPlaneOrientation(int theNum){
172   return myBasePlane[theNum];
173 }
174
175 float VISU_CutPlanesPL::GetRotateX(int theNum){
176   switch(myBasePlane[theNum]){
177   case XY: return myAng[theNum][0];
178   case YZ: return myAng[theNum][1];
179   case ZX: return myAng[theNum][2];
180   }
181   return 0;
182 }
183
184 float VISU_CutPlanesPL::GetRotateY(int theNum){
185   switch(myBasePlane[theNum]){
186   case XY: return myAng[theNum][1];
187   case YZ: return myAng[theNum][2];
188   case ZX: return myAng[theNum][0];
189   }
190   return 0;
191 }
192
193
194 void VISU_CutPlanesPL::SetNbParts(int theNb) {
195   myNbParts = theNb;
196   myPartPosition.resize(myNbParts);
197   myPartCondition.resize(myNbParts,1);
198   Modified();
199 }
200
201
202 void VISU_CutPlanesPL::SetPartPosition(int thePartNumber, float thePartPosition){
203   if(thePartNumber >= myNbParts) return;
204   myPartPosition[thePartNumber] = thePartPosition;
205   myPartCondition[thePartNumber] = 0;
206   Modified();
207 }
208 float VISU_CutPlanesPL::GetPartPosition(int thePartNumber, int theNum){
209   if(thePartNumber >= myNbParts) return 0;
210   float aPosition = myPartPosition[thePartNumber];
211   if(myPartCondition[thePartNumber]){
212       float aDir[3], aBounds[6], aBoundPrj[3];
213       GetInput2()->GetBounds(aBounds);
214       GetDir(aDir,myAng[theNum],myBasePlane[theNum]);
215       GetBoundProject(aBoundPrj,aBounds,aDir);
216       if (myNbParts > 1){
217         float aDBoundPrj = aBoundPrj[2]/(myNbParts - 1);
218         float aDisplacement = aDBoundPrj * myDisplacement[theNum];
219         float aStartPosition = aBoundPrj[0] - 0.5*aDBoundPrj + aDisplacement;
220         aPosition = aStartPosition + thePartNumber*aDBoundPrj;
221       }else
222         aPosition = aBoundPrj[0] + aBoundPrj[2]*myDisplacement[theNum];
223   }
224   return aPosition;
225 }
226
227
228 void VISU_CutPlanesPL::SetPartDefault(int thePartNumber){
229   if(thePartNumber >= myNbParts) return;
230   myPartPosition[thePartNumber] = GetPartPosition(thePartNumber);
231   myPartCondition[thePartNumber] = 1;
232   Modified();
233 }
234 int VISU_CutPlanesPL::IsPartDefault(int thePartNumber){
235   if(thePartNumber >= myNbParts) return 1;
236   return myPartCondition[thePartNumber];
237 }
238
239
240 void VISU_CutPlanesPL::GetDir(float theDir[3],
241                               const float theAng[3],
242                               const PlaneOrientation& theBasePlane)
243 {
244   int iPlane = 0;
245   float aRx[3][3], aRy[3][3], aRz[3][3], aRotation[3][3];
246   switch(theBasePlane){
247   case XY:
248     if(fabs(theAng[0]) > EPS) GetRx(aRx,theAng[0]); else vtkMath::Identity3x3(aRx);
249     if(fabs(theAng[1]) > EPS) GetRy(aRy,theAng[1]); else vtkMath::Identity3x3(aRy);
250     vtkMath::Multiply3x3(aRx,aRy,aRotation);
251     iPlane = 2;
252     break;
253   case YZ:
254     if(fabs(theAng[1]) > EPS) GetRy(aRy,theAng[1]); else vtkMath::Identity3x3(aRy);
255     if(fabs(theAng[2]) > EPS) GetRz(aRz,theAng[2]); else vtkMath::Identity3x3(aRz);
256     vtkMath::Multiply3x3(aRy,aRz,aRotation);
257     iPlane = 0;
258     break;
259   case ZX:
260     if(fabs(theAng[2]) > EPS) GetRz(aRz,theAng[2]); else vtkMath::Identity3x3(aRz);
261     if(fabs(theAng[0]) > EPS) GetRx(aRx,theAng[0]); else vtkMath::Identity3x3(aRx);
262     vtkMath::Multiply3x3(aRz,aRx,aRotation);
263     iPlane = 1;
264     break;
265   }
266   for(int i = 0; i < 3; i++)
267     theDir[i] = aRotation[i][iPlane];
268 }
269
270
271 void VISU_CutPlanesPL::CutWithPlane(vtkAppendPolyData* theAppendPolyData,
272                                     vtkDataSet* theDataSet,
273                                     float theDir[3], float theOrig[3])
274 {
275   vtkCutter *aCutPlane = vtkCutter::New();
276   aCutPlane->SetInput(theDataSet);
277   vtkPlane *aPlane = vtkPlane::New();
278   aPlane->SetOrigin(theOrig);
279
280   aPlane->SetNormal(theDir);
281   aCutPlane->SetCutFunction(aPlane);
282   aPlane->Delete();
283   theAppendPolyData->AddInput(aCutPlane->GetOutput());
284   aCutPlane->Register(theAppendPolyData);
285   aCutPlane->Delete();
286 }
287
288
289 void VISU_CutPlanesPL::CutWithPlanes(vtkAppendPolyData* theAppendPolyData, vtkDataSet* theDataSet,
290                                      int theNbPlanes, float theDir[3], float theBounds[6],
291                                      const vector<float>& thePlanePosition,
292                                      const vector<int>& thePlaneCondition,
293                                      float theDisplacement)
294 {
295   float aBoundPrj[3], aOrig[3], aPosition;
296   GetBoundProject(aBoundPrj, theBounds, theDir);
297   if(theNbPlanes > 1){
298     float aDBoundPrj = aBoundPrj[2]/(theNbPlanes - 1);
299     float aDisplacement = aDBoundPrj*theDisplacement;
300     float aStartPosition = aBoundPrj[0] - 0.5*aDBoundPrj + aDisplacement;
301     for (int i = 0; i < theNbPlanes; i++){
302       aPosition = aStartPosition + i*aDBoundPrj;
303       if(thePlaneCondition[i]){
304         aPosition = aStartPosition + i*aDBoundPrj;
305       }else
306         aPosition = thePlanePosition[i];
307       VISU::Mul(theDir,aPosition,aOrig);
308       CutWithPlane(theAppendPolyData,theDataSet,theDir,aOrig);
309     }
310   }else{
311     if(thePlaneCondition[0])
312       aPosition = aBoundPrj[0] + aBoundPrj[2]*theDisplacement;
313     else
314       aPosition = thePlanePosition[0];
315     VISU::Mul(theDir,aPosition,aOrig);
316     CutWithPlane(theAppendPolyData,theDataSet,theDir,aOrig);
317   }
318   vtkPolyData *aPolyData = theAppendPolyData->GetOutput();
319   aPolyData->Update();
320   theAppendPolyData->Update();
321 }