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