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