Salome HOME
da2d79f9260c10b73b51e9f62888d47edc457b80
[modules/visu.git] / src / PIPELINE / VISU_CutPlanesPL.cxx
1 //  Copyright (C) 2007-2010  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 //  Copyright (C) 2003-2007  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.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22
23 //  VISU OBJECT : interactive object for VISU entities implementation
24 // File:    VISU_PipeLine.cxx
25 // Author:  Alexey PETROV
26 // Module : VISU
27 //
28 #include "VISU_CutPlanesPL.hxx"
29 #include "VISU_FieldTransform.hxx"
30 #include "VISU_PipeLineUtils.hxx"
31 #include "VTKViewer_GeometryFilter.h"
32 #include "VISU_MapperHolder.hxx"
33 #include "VISU_DeformationPL.hxx"
34
35 #include <vtkAppendPolyData.h>
36 #include <vtkCutter.h>
37 #include <vtkPlane.h>
38
39 //#include <vtkUnstructuredGrid.h>
40
41 static vtkFloatingPointType EPS = 1.0E-3;
42
43 #ifdef _DEBUG_
44 static int MYDEBUG = 0;
45 #else
46 static int MYDEBUG = 0;
47 #endif
48
49
50 //----------------------------------------------------------------------------
51 vtkStandardNewMacro(VISU_CutPlanesPL);
52
53
54 //----------------------------------------------------------------------------
55 VISU_CutPlanesPL
56 ::VISU_CutPlanesPL():
57   VISU_OptionalDeformationPL()
58 {
59   if(MYDEBUG) MESSAGE("VISU_CutPlanesPL()::VISU_CutPlanesPL() - "<<this);
60   
61   SetIsShrinkable(false);
62   SetIsFeatureEdgesAllowed(false);
63
64   SetElnoDisassembleState( true );
65
66   myAppendPolyData = vtkAppendPolyData::New();
67
68   myNbParts = 10;
69
70   myBasePlane[0] = YZ;
71   myBasePlane[0] = XY;
72
73   myDisplacement[0] = myDisplacement[1] = 0.5;
74
75   myAng[0][0] = myAng[0][1] = myAng[0][2] = 0.0;
76   myAng[1][0] = myAng[1][1] = myAng[1][2] = 0.0;
77   UseDeformation(false);
78 }
79
80
81 //----------------------------------------------------------------------------
82 VISU_CutPlanesPL
83 ::~VISU_CutPlanesPL()
84 {
85   if(MYDEBUG) MESSAGE("VISU_CutPlanesPL()::~VISU_CutPlanesPL() - "<<this);
86   myAppendPolyData->Delete();
87   myAppendPolyData = NULL;
88 }
89
90
91 //----------------------------------------------------------------------------
92 unsigned long int 
93 VISU_CutPlanesPL
94 ::GetMTime()
95 {
96   unsigned long int aTime = Superclass::GetMTime();
97   
98   if(IsDeformed()) {
99     aTime = std::max(aTime, VISU_OptionalDeformationPL::GetMTime());
100   }
101   
102   aTime = std::max(aTime, myAppendPolyData->GetMTime());
103
104   return aTime;
105 }
106
107
108 //----------------------------------------------------------------------------
109 void
110 VISU_CutPlanesPL
111 ::DoShallowCopy(VISU_PipeLine *thePipeLine,
112                 bool theIsCopyInput)
113 {
114   Superclass::DoShallowCopy(thePipeLine, theIsCopyInput);
115
116   if(VISU_CutPlanesPL *aPipeLine = dynamic_cast<VISU_CutPlanesPL*>(thePipeLine)){
117
118     SetOrientation(aPipeLine->GetPlaneOrientation(),
119                    aPipeLine->GetRotateX(),
120                    aPipeLine->GetRotateY());
121
122     SetDisplacement(aPipeLine->GetDisplacement());
123
124     SetNbParts(aPipeLine->GetNbParts());
125     for (int i = 0, iEnd = GetNbParts(); i < iEnd; i++)
126       if(!aPipeLine->IsPartDefault(i))  
127         SetPartPosition(i, aPipeLine->GetPartPosition(i));
128   }
129 }
130
131
132 //----------------------------------------------------------------------------
133 void
134 VISU_CutPlanesPL
135 ::Init()
136 {
137   Superclass::Init();
138   SetNbParts(10);
139   myBasePlane[0] = YZ;
140   myDisplacement[0] = 0.5;
141   myAng[0][0] = myAng[0][1] = myAng[0][2] = 0.0;
142   SetScale(VISU_DeformationPL::GetDefaultScaleFactor(this));
143 }
144
145
146 //----------------------------------------------------------------------------
147 vtkDataSet* 
148 VISU_CutPlanesPL
149 ::InsertCustomPL()
150 {
151   return GetWarpVectorOutput();
152 }
153
154
155 //----------------------------------------------------------------------------
156 void
157 VISU_CutPlanesPL
158 ::Update()
159 {
160   vtkDataSet* aMergedInput = GetMergedInput();
161   if(VISU::IsQuadraticData(aMergedInput)) // Bug 0020123, note 0005343
162     throw std::runtime_error("Impossible to build presentation");
163
164   ClearAppendPolyData(myAppendPolyData);
165
166
167   if(!myVectorialField || !IsDeformed()){
168     SetMergeFilterInput(aMergedInput,aMergedInput);
169   }
170   
171
172   if(VISU::IsDataOnCells(aMergedInput))
173     GetMapper()->SetScalarModeToUseCellData();
174   else
175     GetMapper()->SetScalarModeToUsePointData();
176
177   SetPartPosition();
178   
179   vtkFloatingPointType aDir[3];
180   GetDir(aDir, 
181          myAng[0], 
182          myBasePlane[0]);
183   
184   vtkFloatingPointType aBounds[6];
185
186   vtkDataSet* aFilterOutput = GetMergeFilterOutput();
187   
188   aFilterOutput->GetBounds(aBounds);
189
190   CutWithPlanes(myAppendPolyData,
191                 aFilterOutput,
192                 myNbParts,
193                 aDir,
194                 aBounds,
195                 myPartPosition, 
196                 myPartCondition, 
197                 myDisplacement[0]);
198
199   
200
201   SetWarpVectorInput(myAppendPolyData->GetOutput());
202   Superclass::Update();
203 }
204
205
206 //----------------------------------------------------------------------------
207 unsigned long int
208 VISU_CutPlanesPL
209 ::GetMemorySize()
210 {
211   unsigned long int aSize = Superclass::GetMemorySize();
212
213   if(vtkDataSet* aDataSet = myAppendPolyData->GetOutput())
214     aSize += aDataSet->GetActualMemorySize() * 1024;
215   
216   int anEnd = myAppendPolyData->GetNumberOfInputConnections(0);
217   for(int anId = 0; anId < anEnd; anId++)
218     if(vtkDataSet* aDataSet = myAppendPolyData->GetInput(anId))
219       aSize += aDataSet->GetActualMemorySize() * 1024;
220
221   return aSize;
222 }
223
224
225 //----------------------------------------------------------------------------
226 void
227 VISU_CutPlanesPL
228 ::SetPartPosition(int theNum)
229 {
230   for(int i = 0; i < myNbParts; i++)
231     myPartPosition[i] = GetPartPosition(i,theNum);
232 }
233
234 void
235 VISU_CutPlanesPL
236 ::ClearAppendPolyData(vtkAppendPolyData *theAppendPolyData)
237 {
238   theAppendPolyData->RemoveAllInputs();
239 }
240
241
242 //----------------------------------------------------------------------------
243 vtkFloatingPointType* 
244 VISU_CutPlanesPL::
245 GetRx(vtkFloatingPointType theRx[3][3], 
246       vtkFloatingPointType thaAng)
247 {
248   theRx[0][0] = 1.0;            theRx[0][1] = 0.0;            theRx[0][2] = 0.0;
249   theRx[1][0] = 0.0;            theRx[1][1] = cos(thaAng);    theRx[1][2] = -sin(thaAng);
250   theRx[2][0] = 0.0;            theRx[2][1] = sin(thaAng);    theRx[2][2] = cos(thaAng);
251   return theRx[0];
252 }
253
254
255
256 //----------------------------------------------------------------------------
257 vtkFloatingPointType* 
258 VISU_CutPlanesPL
259 ::GetRy(vtkFloatingPointType theRy[3][3], 
260         vtkFloatingPointType thaAng)
261 {
262   theRy[0][0] = cos(thaAng);    theRy[0][1] = 0.0;            theRy[0][2] = sin(thaAng);
263   theRy[1][0] = 0.0;            theRy[1][1] = 1.0;            theRy[1][2] = 0.0;
264   theRy[2][0] = -sin(thaAng);   theRy[2][1] = 0.0;            theRy[2][2] = cos(thaAng);
265   return theRy[0];
266 }
267
268
269 //----------------------------------------------------------------------------
270 vtkFloatingPointType* 
271 VISU_CutPlanesPL
272 ::GetRz(vtkFloatingPointType theRz[3][3], 
273         vtkFloatingPointType thaAng)
274 {
275   theRz[0][0] = cos(thaAng);    theRz[0][1] = -sin(thaAng);   theRz[0][2] = 0.0;
276   theRz[1][0] = sin(thaAng);    theRz[1][1] = cos(thaAng);    theRz[1][2] = 0.0;
277   theRz[2][0] = 0.0;            theRz[2][1] = 0.0;            theRz[2][2] = 1.0;
278   return theRz[0];
279 }
280
281
282 //----------------------------------------------------------------------------
283 void
284 VISU_CutPlanesPL
285 ::CorrectPnt(vtkFloatingPointType thePnt[3], 
286              const vtkFloatingPointType BoundPrj[6])
287 {
288   for(int i = 0, j = 0; i < 3; ++i, j=2*i){
289     if(thePnt[i] < BoundPrj[j]) thePnt[i] = BoundPrj[j];
290     if(thePnt[i] > BoundPrj[j+1]) thePnt[i] = BoundPrj[j+1];
291   }
292 }
293
294
295 //----------------------------------------------------------------------------
296 void
297 VISU_CutPlanesPL
298 ::GetBoundProject(vtkFloatingPointType BoundPrj[3], 
299                   const vtkFloatingPointType BoundBox[6], 
300                   const vtkFloatingPointType Dir[3])
301 {
302   vtkFloatingPointType BoundPoints[8][3] = { {BoundBox[0],BoundBox[2],BoundBox[4]},
303                                              {BoundBox[1],BoundBox[2],BoundBox[4]},
304                                              {BoundBox[0],BoundBox[3],BoundBox[4]},
305                                              {BoundBox[1],BoundBox[3],BoundBox[4]},
306                                              {BoundBox[0],BoundBox[2],BoundBox[5]},
307                                              {BoundBox[1],BoundBox[2],BoundBox[5]},
308                                              {BoundBox[0],BoundBox[3],BoundBox[5]},
309                                              {BoundBox[1],BoundBox[3],BoundBox[5]}};
310   BoundPrj[0] = vtkMath::Dot(Dir,BoundPoints[0]), BoundPrj[1] = BoundPrj[0];
311   for(int i = 1; i < 8; i++){
312     vtkFloatingPointType tmp = vtkMath::Dot(Dir,BoundPoints[i]);
313     if(BoundPrj[1] < tmp) BoundPrj[1] = tmp;
314     if(BoundPrj[0] > tmp) BoundPrj[0] = tmp;
315   }
316   BoundPrj[2] = BoundPrj[1] - BoundPrj[0];
317   BoundPrj[1] = BoundPrj[0] + (1.0 - EPS)*BoundPrj[2];
318   BoundPrj[0] = BoundPrj[0] + EPS*BoundPrj[2];
319   BoundPrj[2] = BoundPrj[1] - BoundPrj[0];
320 }
321
322
323 //----------------------------------------------------------------------------
324 void
325 VISU_CutPlanesPL
326 ::SetOrientation(const VISU_CutPlanesPL::PlaneOrientation& theOrient,
327                  vtkFloatingPointType theXAng, 
328                  vtkFloatingPointType theYAng, 
329                  int theNum)
330 {
331   myBasePlane[theNum] = theOrient;
332   switch(myBasePlane[theNum]){
333   case XY: myAng[theNum][0] = theXAng; break;
334   case YZ: myAng[theNum][1] = theXAng; break;
335   case ZX: myAng[theNum][2] = theXAng; break;
336   }
337   switch(myBasePlane[theNum]){
338   case XY: myAng[theNum][1] = theYAng; break;
339   case YZ: myAng[theNum][2] = theYAng; break;
340   case ZX: myAng[theNum][0] = theYAng; break;
341   }
342 }
343
344
345 //----------------------------------------------------------------------------
346 const VISU_CutPlanesPL::PlaneOrientation& 
347 VISU_CutPlanesPL
348 ::GetPlaneOrientation(int theNum)
349 {
350   return myBasePlane[theNum];
351 }
352
353 vtkFloatingPointType 
354 VISU_CutPlanesPL
355 ::GetRotateX(int theNum)
356 {
357   switch(myBasePlane[theNum]){
358   case XY: return myAng[theNum][0];
359   case YZ: return myAng[theNum][1];
360   case ZX: return myAng[theNum][2];
361   }
362   return 0;
363 }
364
365
366 //----------------------------------------------------------------------------
367 vtkFloatingPointType 
368 VISU_CutPlanesPL
369 ::GetRotateY(int theNum)
370 {
371   switch(myBasePlane[theNum]){
372   case XY: return myAng[theNum][1];
373   case YZ: return myAng[theNum][2];
374   case ZX: return myAng[theNum][0];
375   }
376   return 0;
377 }
378
379
380 //----------------------------------------------------------------------------
381 vtkFloatingPointType 
382 VISU_CutPlanesPL
383 ::GetDisplacement(int theNum)
384 {
385   return myDisplacement[theNum];
386 }
387
388
389 //----------------------------------------------------------------------------
390 void
391 VISU_CutPlanesPL
392 ::SetDisplacement(vtkFloatingPointType theDisp, 
393                   int theNum) 
394 {
395   if(VISU::CheckIsSameValue(myDisplacement[theNum], theDisp))
396     return;
397
398   myDisplacement[theNum] = theDisp;
399   Modified();
400 }
401
402
403 //----------------------------------------------------------------------------
404 void
405 VISU_CutPlanesPL
406 ::SetNbParts(int theNbParts) 
407 {
408   if(theNbParts > 0 && GetNbParts() != theNbParts){
409     myPartPosition.resize(theNbParts);
410     myPartCondition.resize(theNbParts, 1);
411     myNbParts = theNbParts;
412     Modified();
413   }
414 }
415
416
417 //----------------------------------------------------------------------------
418 int
419 VISU_CutPlanesPL
420 ::GetNbParts() 
421 {
422   return myPartPosition.size();
423 }
424
425
426 //----------------------------------------------------------------------------
427 void
428 VISU_CutPlanesPL
429 ::SetPartPosition(int thePartNumber, 
430                   vtkFloatingPointType thePartPosition)
431 {
432   if(thePartNumber >= myNbParts) 
433     return;
434
435   bool anIsSameValue = VISU::CheckIsSameValue(myPartPosition[thePartNumber], thePartPosition);
436   anIsSameValue &= VISU::CheckIsSameValue(myPartCondition[thePartNumber], 0);
437   if(anIsSameValue)
438     return;
439
440   myPartPosition[thePartNumber] = thePartPosition;
441   myPartCondition[thePartNumber] = 0;
442   Modified();
443 }
444
445
446 //----------------------------------------------------------------------------
447 vtkFloatingPointType 
448 VISU_CutPlanesPL
449 ::GetPartPosition(int thePartNumber, 
450                   int theNum)
451 {
452   if(thePartNumber >= myNbParts) 
453     return 0;
454
455   vtkFloatingPointType aPosition = myPartPosition[thePartNumber];
456   if(myPartCondition[thePartNumber]){
457       vtkFloatingPointType aDir[3], aBounds[6], aBoundPrj[3];
458       if(!IsDeformed()) 
459         GetMergedInput()->GetBounds(aBounds);
460       else
461         GetMergeFilterOutput()->GetBounds(aBounds);
462
463
464       GetDir(aDir,
465              myAng[theNum],
466              myBasePlane[theNum]);
467
468       GetBoundProject(aBoundPrj,
469                       aBounds,
470                       aDir);
471
472       if(myNbParts > 1){
473         vtkFloatingPointType aDBoundPrj = aBoundPrj[2]/(myNbParts - 1);
474         vtkFloatingPointType aDisplacement = aDBoundPrj * myDisplacement[theNum];
475         vtkFloatingPointType aStartPosition = aBoundPrj[0] - 0.5*aDBoundPrj + aDisplacement;
476         aPosition = aStartPosition + thePartNumber*aDBoundPrj;
477       }else
478         aPosition = aBoundPrj[0] + aBoundPrj[2]*myDisplacement[theNum];
479   }
480
481   return aPosition;
482 }
483
484
485 //----------------------------------------------------------------------------
486 void
487 VISU_CutPlanesPL
488 ::SetPartDefault(int thePartNumber)
489 {
490   if(thePartNumber >= myNbParts) 
491     return;
492
493   bool anIsSameValue = VISU::CheckIsSameValue(myPartPosition[thePartNumber], GetPartPosition(thePartNumber));
494   anIsSameValue &= VISU::CheckIsSameValue(myPartCondition[thePartNumber], 1);
495   if(anIsSameValue)
496     return;
497
498   myPartPosition[thePartNumber] = GetPartPosition(thePartNumber);
499   myPartCondition[thePartNumber] = 1;
500   Modified();
501 }
502
503
504 //----------------------------------------------------------------------------
505 int
506 VISU_CutPlanesPL
507 ::IsPartDefault(int thePartNumber)
508 {
509   if(thePartNumber >= myNbParts) 
510     return 1;
511
512   return myPartCondition[thePartNumber];
513 }
514
515
516 //----------------------------------------------------------------------------
517 void
518 VISU_CutPlanesPL
519 ::GetDir(vtkFloatingPointType theDir[3],
520          const vtkFloatingPointType theAng[3],
521          const PlaneOrientation& theBasePlane)
522 {
523   int iPlane = 0;
524   vtkFloatingPointType aRx[3][3], aRy[3][3], aRz[3][3], aRotation[3][3];
525   switch(theBasePlane){
526   case XY:
527     if(fabs(theAng[0]) > EPS) GetRx(aRx,theAng[0]); else vtkMath::Identity3x3(aRx);
528     if(fabs(theAng[1]) > EPS) GetRy(aRy,theAng[1]); else vtkMath::Identity3x3(aRy);
529     vtkMath::Multiply3x3(aRx,aRy,aRotation);
530     iPlane = 2;
531     break;
532   case YZ:
533     if(fabs(theAng[1]) > EPS) GetRy(aRy,theAng[1]); else vtkMath::Identity3x3(aRy);
534     if(fabs(theAng[2]) > EPS) GetRz(aRz,theAng[2]); else vtkMath::Identity3x3(aRz);
535     vtkMath::Multiply3x3(aRy,aRz,aRotation);
536     iPlane = 0;
537     break;
538   case ZX:
539     if(fabs(theAng[2]) > EPS) GetRz(aRz,theAng[2]); else vtkMath::Identity3x3(aRz);
540     if(fabs(theAng[0]) > EPS) GetRx(aRx,theAng[0]); else vtkMath::Identity3x3(aRx);
541     vtkMath::Multiply3x3(aRz,aRx,aRotation);
542     iPlane = 1;
543     break;
544   }
545
546   for(int i = 0; i < 3; i++)
547     theDir[i] = aRotation[i][iPlane];
548 }
549
550
551 //----------------------------------------------------------------------------
552 void
553 VISU_CutPlanesPL
554 ::CutWithPlane(vtkAppendPolyData* theAppendPolyData,
555                vtkDataSet* theDataSet,
556                vtkFloatingPointType theDir[3], 
557                vtkFloatingPointType theOrig[3])
558 {
559   vtkCutter *aCutPlane = vtkCutter::New();
560   aCutPlane->SetInput(theDataSet);
561   vtkPlane *aPlane = vtkPlane::New();
562   aPlane->SetOrigin(theOrig);
563
564   aPlane->SetNormal(theDir);
565   aCutPlane->SetCutFunction(aPlane);
566   aPlane->Delete();
567   theAppendPolyData->AddInput(aCutPlane->GetOutput());
568   aCutPlane->Delete();
569 }
570
571
572 //----------------------------------------------------------------------------
573 void
574 VISU_CutPlanesPL
575 ::CutWithPlanes(vtkAppendPolyData* theAppendPolyData, 
576                 vtkDataSet* theDataSet,
577                 int theNbPlanes, 
578                 vtkFloatingPointType theDir[3], 
579                 vtkFloatingPointType theBounds[6],
580                 const std::vector<vtkFloatingPointType>& thePlanePosition,
581                 const std::vector<int>& thePlaneCondition,
582                 vtkFloatingPointType theDisplacement)
583 {
584   vtkFloatingPointType aBoundPrj[3], aOrig[3], aPosition;
585   GetBoundProject(aBoundPrj, theBounds, theDir);
586   if(theNbPlanes > 1){
587     vtkFloatingPointType aDBoundPrj = aBoundPrj[2]/(theNbPlanes - 1);
588     vtkFloatingPointType aDisplacement = aDBoundPrj*theDisplacement;
589     vtkFloatingPointType aStartPosition = aBoundPrj[0] - 0.5*aDBoundPrj + aDisplacement;
590     for (int i = 0; i < theNbPlanes; i++){
591       aPosition = aStartPosition + i*aDBoundPrj;
592       if(thePlaneCondition[i]){
593         aPosition = aStartPosition + i*aDBoundPrj;
594       }else
595         aPosition = thePlanePosition[i];
596       VISU::Mul(theDir,aPosition,aOrig);
597       CutWithPlane(theAppendPolyData,theDataSet,theDir,aOrig);
598     }
599   }else{
600     if(thePlaneCondition[0])
601       aPosition = aBoundPrj[0] + aBoundPrj[2]*theDisplacement;
602     else
603       aPosition = thePlanePosition[0];
604     VISU::Mul(theDir,aPosition,aOrig);
605     CutWithPlane(theAppendPolyData,theDataSet,theDir,aOrig);
606   }
607   vtkPolyData *aPolyData = theAppendPolyData->GetOutput();
608   aPolyData->Update();
609   theAppendPolyData->Update();
610 }
611
612
613 //----------------------------------------------------------------------------
614 void
615 VISU_CutPlanesPL::SetVectorialField(VISU::PUnstructuredGridIDMapper theMapper)
616 {  
617   if(myVectorialField == theMapper)
618     return;
619
620   if(CheckCanDeformate(theMapper->GetOutput())){
621     myVectorialField = theMapper;
622     
623     SetMergeFilterInput(GetMergedInput(),theMapper->GetOutput());
624   }
625   else
626     UseDeformation(false);
627   
628   Modified();
629 }
630
631 //----------------------------------------------------------------------------
632 VISU::PUnstructuredGridIDMapper VISU_CutPlanesPL::
633 getVectorialField()
634 {
635   return myVectorialField;
636 }
637
638 //----------------------------------------------------------------------------
639 void VISU_CutPlanesPL::SetMapScale(vtkFloatingPointType theMapScale){
640   Superclass::SetMapScale(theMapScale);
641   if(IsDeformed())
642     VISU_OptionalDeformationPL::SetMapScale(theMapScale);
643 }