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 SetPlaneOrientation(aPipeLine->GetPlaneOrientation());
54 SetDisplacement(aPipeLine->GetDisplacement());
55 SetNbParts(aPipeLine->GetNbParts());
56 for (int i = 0, iend = GetNbParts(); i < iend; i++)
57 if(!aPipeLine->IsPartDefault(i)) SetPartPosition(i, aPipeLine->GetPartPosition(i));
58 SetRotateX(aPipeLine->GetRotateX());
59 SetRotateY(aPipeLine->GetRotateY());
63 void VISU_CutPlanesPL::Init(){
64 VISU_ScalarMapPL::Init();
68 myDisplacement[0] = 0.5;
69 myAng[0][0] = myAng[0][1] = myAng[0][2] = 0.0;
72 VISU_ScalarMapPL::THook* VISU_CutPlanesPL::DoHook(){
73 return myAppendPolyData->GetOutput();
76 void VISU_CutPlanesPL::Update(){
77 ClearAppendPolyData(myAppendPolyData);
80 GetDir(aDir,myAng[0],myBasePlane[0]);
82 myInput->GetBounds(aBounds);
83 vtkDataSet* aDataSet = myFieldTransform->GetUnstructuredGridOutput();
84 CutWithPlanes(myAppendPolyData,aDataSet,myNbParts,aDir,aBounds,
85 myPartPosition,myPartCondition,myDisplacement[0]);
87 VISU_ScalarMapPL::Update();
90 void VISU_CutPlanesPL::SetPartPosition(int theNum){
91 for(int i = 0; i < myNbParts; i++)
92 myPartPosition[i] = GetPartPosition(i,theNum);
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));
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);
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);
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;
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];
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;
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];
153 void VISU_CutPlanesPL::SetRotateX(float theAng, int theNum){
154 switch(myBasePlane[theNum]){
155 case XY: myAng[theNum][0] = theAng; break;
156 case YZ: myAng[theNum][1] = theAng; break;
157 case ZX: myAng[theNum][2] = theAng; break;
161 float VISU_CutPlanesPL::GetRotateX(int theNum){
162 switch(myBasePlane[theNum]){
163 case XY: return myAng[theNum][0];
164 case YZ: return myAng[theNum][1];
165 case ZX: return myAng[theNum][2];
169 void VISU_CutPlanesPL::SetRotateY(float theAng, int theNum){
170 switch(myBasePlane[theNum]){
171 case XY: myAng[theNum][1] = theAng; break;
172 case YZ: myAng[theNum][2] = theAng; break;
173 case ZX: myAng[theNum][0] = theAng; break;
177 float VISU_CutPlanesPL::GetRotateY(int theNum){
178 switch(myBasePlane[theNum]){
179 case XY: return myAng[theNum][1];
180 case YZ: return myAng[theNum][2];
181 case ZX: return myAng[theNum][0];
186 void VISU_CutPlanesPL::SetNbParts(int theNb) {
188 myPartPosition.resize(myNbParts);
189 myPartCondition.resize(myNbParts,1);
194 void VISU_CutPlanesPL::SetPartPosition(int thePartNumber, float thePartPosition){
195 if(thePartNumber >= myNbParts) return;
196 myPartPosition[thePartNumber] = thePartPosition;
197 myPartCondition[thePartNumber] = 0;
200 float VISU_CutPlanesPL::GetPartPosition(int thePartNumber, int theNum){
201 if(thePartNumber >= myNbParts) return 0;
202 float aPosition = myPartPosition[thePartNumber];
203 if(myPartCondition[thePartNumber]){
204 float aDir[3], aBounds[6], aBoundPrj[3];
205 myInput->GetBounds(aBounds);
206 GetDir(aDir,myAng[theNum],myBasePlane[theNum]);
207 GetBoundProject(aBoundPrj,aBounds,aDir);
209 float aDBoundPrj = aBoundPrj[2]/(myNbParts - 1);
210 float aDisplacement = aDBoundPrj * myDisplacement[theNum];
211 float aStartPosition = aBoundPrj[0] - 0.5*aDBoundPrj + aDisplacement;
212 aPosition = aStartPosition + thePartNumber*aDBoundPrj;
214 aPosition = aBoundPrj[0] + aBoundPrj[2]*myDisplacement[theNum];
220 void VISU_CutPlanesPL::SetPartDefault(int thePartNumber){
221 if(thePartNumber >= myNbParts) return;
222 myPartPosition[thePartNumber] = GetPartPosition(thePartNumber);
223 myPartCondition[thePartNumber] = 1;
226 int VISU_CutPlanesPL::IsPartDefault(int thePartNumber){
227 if(thePartNumber >= myNbParts) return 1;
228 return myPartCondition[thePartNumber];
232 void VISU_CutPlanesPL::GetDir(float theDir[3], float theAng[3], const PlaneOrientation& theBasePlane){
234 float aRx[3][3], aRy[3][3], aRz[3][3], aRotation[3][3];
235 switch(theBasePlane){
237 if(fabs(theAng[0]) > EPS) GetRx(aRx,theAng[0]); else vtkMath::Identity3x3(aRx);
238 if(fabs(theAng[1]) > EPS) GetRy(aRy,theAng[1]); else vtkMath::Identity3x3(aRy);
239 vtkMath::Multiply3x3(aRx,aRy,aRotation);
243 if(fabs(theAng[1]) > EPS) GetRy(aRy,theAng[1]); else vtkMath::Identity3x3(aRy);
244 if(fabs(theAng[2]) > EPS) GetRz(aRz,theAng[2]); else vtkMath::Identity3x3(aRz);
245 vtkMath::Multiply3x3(aRy,aRz,aRotation);
249 if(fabs(theAng[2]) > EPS) GetRz(aRz,theAng[2]); else vtkMath::Identity3x3(aRz);
250 if(fabs(theAng[0]) > EPS) GetRx(aRx,theAng[0]); else vtkMath::Identity3x3(aRx);
251 vtkMath::Multiply3x3(aRz,aRx,aRotation);
255 for(int i = 0; i < 3; i++)
256 theDir[i] = aRotation[i][iPlane];
260 void VISU_CutPlanesPL::CutWithPlane(vtkAppendPolyData* theAppendPolyData,
261 vtkDataSet* theDataSet,
262 float theDir[3], float theOrig[3])
264 vtkCutter *aCutPlane = vtkCutter::New();
265 aCutPlane->SetInput(theDataSet);
266 vtkPlane *aPlane = vtkPlane::New();
267 aPlane->SetOrigin(theOrig);
269 aPlane->SetNormal(theDir);
270 aCutPlane->SetCutFunction(aPlane);
272 theAppendPolyData->AddInput(aCutPlane->GetOutput());
273 aCutPlane->Register(theAppendPolyData);
278 void VISU_CutPlanesPL::CutWithPlanes(vtkAppendPolyData* theAppendPolyData, vtkDataSet* theDataSet,
279 int theNbPlanes, float theDir[3], float theBounds[6],
280 const vector<float>& thePlanePosition,
281 const vector<int>& thePlaneCondition,
282 float theDisplacement)
284 float aBoundPrj[3], aOrig[3], aPosition;
285 GetBoundProject(aBoundPrj, theBounds, theDir);
287 float aDBoundPrj = aBoundPrj[2]/(theNbPlanes - 1);
288 float aDisplacement = aDBoundPrj*theDisplacement;
289 float aStartPosition = aBoundPrj[0] - 0.5*aDBoundPrj + aDisplacement;
290 for (int i = 0; i < theNbPlanes; i++){
291 aPosition = aStartPosition + i*aDBoundPrj;
292 float aDelta = (aBoundPrj[0] - aPosition) / aBoundPrj[2];
293 if(thePlaneCondition[i]){
294 aPosition = aStartPosition + i*aDBoundPrj;
296 aPosition = thePlanePosition[i];
297 VISU::Mul(theDir,aPosition,aOrig);
298 CutWithPlane(theAppendPolyData,theDataSet,theDir,aOrig);
301 if(thePlaneCondition[0])
302 aPosition = aBoundPrj[0] + aBoundPrj[2]*theDisplacement;
304 aPosition = thePlanePosition[0];
305 VISU::Mul(theDir,aPosition,aOrig);
306 CutWithPlane(theAppendPolyData,theDataSet,theDir,aOrig);
308 vtkPolyData *aPolyData = theAppendPolyData->GetOutput();
310 theAppendPolyData->Update();