1 // VISU OBJECT : interactive object for VISU entities implementation
3 // Copyright (C) 2003 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
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.
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.
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
20 // See http://www.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org
23 // File: VISU_PipeLine.cxx
24 // Author: Alexey PETROV
28 #include "VISU_CutPlanesPL.hxx"
29 #include "VISU_PipeLineUtils.hxx"
30 #include "SALOME_GeometryFilter.h"
32 #include <vtkAppendPolyData.h>
33 #include <vtkCutter.h>
38 static float EPS = 1.0E-3;
40 vtkStandardNewMacro(VISU_CutPlanesPL);
42 VISU_CutPlanesPL::VISU_CutPlanesPL(){
43 myAppendPolyData = vtkAppendPolyData::New();
46 VISU_CutPlanesPL::~VISU_CutPlanesPL(){
47 myAppendPolyData->Delete();
50 void VISU_CutPlanesPL::ShallowCopy(VISU_PipeLine *thePipeLine){
51 VISU_ScalarMapPL::ShallowCopy(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));
62 void VISU_CutPlanesPL::Init(){
63 VISU_ScalarMapPL::Init();
67 myDisplacement[0] = 0.5;
68 myAng[0][0] = myAng[0][1] = myAng[0][2] = 0.0;
71 VISU_ScalarMapPL::THook* VISU_CutPlanesPL::DoHook(){
72 return myAppendPolyData->GetOutput();
75 void VISU_CutPlanesPL::Update(){
76 ClearAppendPolyData(myAppendPolyData);
79 GetDir(aDir,myAng[0],myBasePlane[0]);
81 myInput->GetBounds(aBounds);
82 vtkDataSet* aDataSet = myFieldTransform->GetUnstructuredGridOutput();
83 CutWithPlanes(myAppendPolyData,aDataSet,myNbParts,aDir,aBounds,
84 myPartPosition,myPartCondition,myDisplacement[0]);
86 VISU_ScalarMapPL::Update();
89 void VISU_CutPlanesPL::SetPartPosition(int theNum){
90 for(int i = 0; i < myNbParts; i++)
91 myPartPosition[i] = GetPartPosition(i,theNum);
94 void VISU_CutPlanesPL::ClearAppendPolyData(vtkAppendPolyData *theAppendPolyData){
95 int iEnd = theAppendPolyData->GetNumberOfInputs();
96 for(int i = iEnd-1; i >= 0; i--)
97 theAppendPolyData->RemoveInput(theAppendPolyData->GetInput(i));
100 float* VISU_CutPlanesPL::GetRx(float theRx[3][3], float thaAng){
101 theRx[0][0] = 1.0; theRx[0][1] = 0.0; theRx[0][2] = 0.0;
102 theRx[1][0] = 0.0; theRx[1][1] = cos(thaAng); theRx[1][2] = -sin(thaAng);
103 theRx[2][0] = 0.0; theRx[2][1] = sin(thaAng); theRx[2][2] = cos(thaAng);
108 float* VISU_CutPlanesPL::GetRy(float theRy[3][3], float thaAng){
109 theRy[0][0] = cos(thaAng); theRy[0][1] = 0.0; theRy[0][2] = sin(thaAng);
110 theRy[1][0] = 0.0; theRy[1][1] = 1.0; theRy[1][2] = 0.0;
111 theRy[2][0] = -sin(thaAng); theRy[2][1] = 0.0; theRy[2][2] = cos(thaAng);
116 float* VISU_CutPlanesPL::GetRz(float theRz[3][3], float thaAng){
117 theRz[0][0] = cos(thaAng); theRz[0][1] = -sin(thaAng); theRz[0][2] = 0.0;
118 theRz[1][0] = sin(thaAng); theRz[1][1] = cos(thaAng); theRz[1][2] = 0.0;
119 theRz[2][0] = 0.0; theRz[2][1] = 0.0; theRz[2][2] = 1.0;
124 void VISU_CutPlanesPL::CorrectPnt(float thePnt[3], const float BoundPrj[6]){
125 for(int i = 0, j = 0; i < 3; ++i, j=2*i){
126 if(thePnt[i] < BoundPrj[j]) thePnt[i] = BoundPrj[j];
127 if(thePnt[i] > BoundPrj[j+1]) thePnt[i] = BoundPrj[j+1];
131 void VISU_CutPlanesPL::GetBoundProject(float BoundPrj[3], const float BoundBox[6], const float Dir[3]){
132 float BoundPoints[8][3] = { {BoundBox[0],BoundBox[2],BoundBox[4]},
133 {BoundBox[1],BoundBox[2],BoundBox[4]},
134 {BoundBox[0],BoundBox[3],BoundBox[4]},
135 {BoundBox[1],BoundBox[3],BoundBox[4]},
136 {BoundBox[0],BoundBox[2],BoundBox[5]},
137 {BoundBox[1],BoundBox[2],BoundBox[5]},
138 {BoundBox[0],BoundBox[3],BoundBox[5]},
139 {BoundBox[1],BoundBox[3],BoundBox[5]}};
140 BoundPrj[0] = vtkMath::Dot(Dir,BoundPoints[0]), BoundPrj[1] = BoundPrj[0];
141 for(int i = 1; i < 8; i++){
142 float tmp = vtkMath::Dot(Dir,BoundPoints[i]);
143 if(BoundPrj[1] < tmp) BoundPrj[1] = tmp;
144 if(BoundPrj[0] > tmp) BoundPrj[0] = tmp;
146 BoundPrj[2] = BoundPrj[1] - BoundPrj[0];
147 BoundPrj[1] = BoundPrj[0] + (1.0 - EPS)*BoundPrj[2];
148 BoundPrj[0] = BoundPrj[0] + EPS*BoundPrj[2];
149 BoundPrj[2] = BoundPrj[1] - BoundPrj[0];
153 void VISU_CutPlanesPL::SetOrientation(const VISU_CutPlanesPL::PlaneOrientation& theOrient,
154 float theXAng, float theYAng, int theNum)
156 myBasePlane[theNum] = theOrient;
157 switch(myBasePlane[theNum]){
158 case XY: myAng[theNum][0] = theXAng; break;
159 case YZ: myAng[theNum][1] = theXAng; break;
160 case ZX: myAng[theNum][2] = theXAng; break;
162 switch(myBasePlane[theNum]){
163 case XY: myAng[theNum][1] = theYAng; break;
164 case YZ: myAng[theNum][2] = theYAng; break;
165 case ZX: myAng[theNum][0] = theYAng; break;
170 const VISU_CutPlanesPL::PlaneOrientation& VISU_CutPlanesPL::GetPlaneOrientation(int theNum){
171 return myBasePlane[theNum];
174 float VISU_CutPlanesPL::GetRotateX(int theNum){
175 switch(myBasePlane[theNum]){
176 case XY: return myAng[theNum][0];
177 case YZ: return myAng[theNum][1];
178 case ZX: return myAng[theNum][2];
182 float VISU_CutPlanesPL::GetRotateY(int theNum){
183 switch(myBasePlane[theNum]){
184 case XY: return myAng[theNum][1];
185 case YZ: return myAng[theNum][2];
186 case ZX: return myAng[theNum][0];
191 void VISU_CutPlanesPL::SetNbParts(int theNb) {
193 myPartPosition.resize(myNbParts);
194 myPartCondition.resize(myNbParts,1);
199 void VISU_CutPlanesPL::SetPartPosition(int thePartNumber, float thePartPosition){
200 if(thePartNumber >= myNbParts) return;
201 myPartPosition[thePartNumber] = thePartPosition;
202 myPartCondition[thePartNumber] = 0;
205 float VISU_CutPlanesPL::GetPartPosition(int thePartNumber, int theNum){
206 if(thePartNumber >= myNbParts) return 0;
207 float aPosition = myPartPosition[thePartNumber];
208 if(myPartCondition[thePartNumber]){
209 float aDir[3], aBounds[6], aBoundPrj[3];
210 myInput->GetBounds(aBounds);
211 GetDir(aDir,myAng[theNum],myBasePlane[theNum]);
212 GetBoundProject(aBoundPrj,aBounds,aDir);
214 float aDBoundPrj = aBoundPrj[2]/(myNbParts - 1);
215 float aDisplacement = aDBoundPrj * myDisplacement[theNum];
216 float aStartPosition = aBoundPrj[0] - 0.5*aDBoundPrj + aDisplacement;
217 aPosition = aStartPosition + thePartNumber*aDBoundPrj;
219 aPosition = aBoundPrj[0] + aBoundPrj[2]*myDisplacement[theNum];
225 void VISU_CutPlanesPL::SetPartDefault(int thePartNumber){
226 if(thePartNumber >= myNbParts) return;
227 myPartPosition[thePartNumber] = GetPartPosition(thePartNumber);
228 myPartCondition[thePartNumber] = 1;
231 int VISU_CutPlanesPL::IsPartDefault(int thePartNumber){
232 if(thePartNumber >= myNbParts) return 1;
233 return myPartCondition[thePartNumber];
237 void VISU_CutPlanesPL::GetDir(float theDir[3], float theAng[3], const PlaneOrientation& theBasePlane){
239 float aRx[3][3], aRy[3][3], aRz[3][3], aRotation[3][3];
240 switch(theBasePlane){
242 if(fabs(theAng[0]) > EPS) GetRx(aRx,theAng[0]); else vtkMath::Identity3x3(aRx);
243 if(fabs(theAng[1]) > EPS) GetRy(aRy,theAng[1]); else vtkMath::Identity3x3(aRy);
244 vtkMath::Multiply3x3(aRx,aRy,aRotation);
248 if(fabs(theAng[1]) > EPS) GetRy(aRy,theAng[1]); else vtkMath::Identity3x3(aRy);
249 if(fabs(theAng[2]) > EPS) GetRz(aRz,theAng[2]); else vtkMath::Identity3x3(aRz);
250 vtkMath::Multiply3x3(aRy,aRz,aRotation);
254 if(fabs(theAng[2]) > EPS) GetRz(aRz,theAng[2]); else vtkMath::Identity3x3(aRz);
255 if(fabs(theAng[0]) > EPS) GetRx(aRx,theAng[0]); else vtkMath::Identity3x3(aRx);
256 vtkMath::Multiply3x3(aRz,aRx,aRotation);
260 for(int i = 0; i < 3; i++)
261 theDir[i] = aRotation[i][iPlane];
265 void VISU_CutPlanesPL::CutWithPlane(vtkAppendPolyData* theAppendPolyData,
266 vtkDataSet* theDataSet,
267 float theDir[3], float theOrig[3])
269 vtkCutter *aCutPlane = vtkCutter::New();
270 aCutPlane->SetInput(theDataSet);
271 vtkPlane *aPlane = vtkPlane::New();
272 aPlane->SetOrigin(theOrig);
274 aPlane->SetNormal(theDir);
275 aCutPlane->SetCutFunction(aPlane);
277 theAppendPolyData->AddInput(aCutPlane->GetOutput());
278 aCutPlane->Register(theAppendPolyData);
283 void VISU_CutPlanesPL::CutWithPlanes(vtkAppendPolyData* theAppendPolyData, vtkDataSet* theDataSet,
284 int theNbPlanes, float theDir[3], float theBounds[6],
285 const vector<float>& thePlanePosition,
286 const vector<int>& thePlaneCondition,
287 float theDisplacement)
289 float aBoundPrj[3], aOrig[3], aPosition;
290 GetBoundProject(aBoundPrj, theBounds, theDir);
292 float aDBoundPrj = aBoundPrj[2]/(theNbPlanes - 1);
293 float aDisplacement = aDBoundPrj*theDisplacement;
294 float aStartPosition = aBoundPrj[0] - 0.5*aDBoundPrj + aDisplacement;
295 for (int i = 0; i < theNbPlanes; i++){
296 aPosition = aStartPosition + i*aDBoundPrj;
297 float aDelta = (aBoundPrj[0] - aPosition) / aBoundPrj[2];
298 if(thePlaneCondition[i]){
299 aPosition = aStartPosition + i*aDBoundPrj;
301 aPosition = thePlanePosition[i];
302 VISU::Mul(theDir,aPosition,aOrig);
303 CutWithPlane(theAppendPolyData,theDataSet,theDir,aOrig);
306 if(thePlaneCondition[0])
307 aPosition = aBoundPrj[0] + aBoundPrj[2]*theDisplacement;
309 aPosition = thePlanePosition[0];
310 VISU::Mul(theDir,aPosition,aOrig);
311 CutWithPlane(theAppendPolyData,theDataSet,theDir,aOrig);
313 vtkPolyData *aPolyData = theAppendPolyData->GetOutput();
315 theAppendPolyData->Update();