]> SALOME platform Git repositories - modules/gui.git/blob - src/VTKViewer/VTKViewer_Utilities.cxx
Salome HOME
6f5a701ff0adc8290ee51e9f9138cfb7cfa8429a
[modules/gui.git] / src / VTKViewer / VTKViewer_Utilities.cxx
1 // Copyright (C) 2005  OPEN CASCADE, CEA/DEN, EDF R&D, PRINCIPIA R&D
2 // 
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either 
6 // version 2.1 of the License.
7 // 
8 // This library is distributed in the hope that it will be useful 
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of 
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public  
14 // License along with this library; if not, write to the Free Software 
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/
18 //
19
20 #include "VTKViewer_Utilities.h"
21 #include "VTKViewer_Actor.h"
22
23 #include <algorithm>
24
25 // VTK Includes
26 #include <vtkMath.h>
27 #include <vtkCamera.h>
28 #include <vtkRenderer.h>
29 #include <vtkRenderWindow.h>
30
31 using namespace std;
32
33 /*!@see vtkRenderer::ResetCamera(vtkFloatingPointType bounds[6]) method*/
34 void 
35 ResetCamera(vtkRenderer* theRenderer, 
36             int theUsingZeroFocalPoint)
37 {  
38   if(!theRenderer)
39     return;
40
41   vtkCamera* aCamera = theRenderer->GetActiveCamera();
42   if(!aCamera) 
43     return;
44
45   vtkFloatingPointType aBounds[6];
46   int aCount = ComputeVisiblePropBounds(theRenderer,aBounds);
47
48   if(theUsingZeroFocalPoint || aCount){
49     static vtkFloatingPointType MIN_DISTANCE = 1.0 / VTK_LARGE_FLOAT;
50
51     vtkFloatingPointType aLength = aBounds[1]-aBounds[0];
52     aLength = max((aBounds[3]-aBounds[2]),aLength);
53     aLength = max((aBounds[5]-aBounds[4]),aLength);
54     
55     if(aLength < MIN_DISTANCE)
56       return;
57
58     vtkFloatingPointType aWidth = 
59       sqrt((aBounds[1]-aBounds[0])*(aBounds[1]-aBounds[0]) +
60            (aBounds[3]-aBounds[2])*(aBounds[3]-aBounds[2]) +
61            (aBounds[5]-aBounds[4])*(aBounds[5]-aBounds[4]));
62     
63     if(aWidth < MIN_DISTANCE)
64       return;
65
66     vtkFloatingPointType aViewPlaneNormal[3];
67     aCamera->GetViewPlaneNormal(aViewPlaneNormal);
68     
69     vtkFloatingPointType aCenter[3] = {0.0, 0.0, 0.0};
70     if(!theUsingZeroFocalPoint){
71       aCenter[0] = (aBounds[0] + aBounds[1])/2.0;
72       aCenter[1] = (aBounds[2] + aBounds[3])/2.0;
73       aCenter[2] = (aBounds[4] + aBounds[5])/2.0;
74     }
75     aCamera->SetFocalPoint(aCenter[0],aCenter[1],aCenter[2]);
76     
77     vtkFloatingPointType aViewAngle = aCamera->GetViewAngle();
78     vtkFloatingPointType aDistance = 2.0*aWidth/tan(aViewAngle*vtkMath::Pi()/360.0);
79     
80     // check view-up vector against view plane normal
81     vtkFloatingPointType aViewUp[3];
82     aCamera->GetViewUp(aViewUp);
83     if(fabs(vtkMath::Dot(aViewUp,aViewPlaneNormal)) > 0.999)
84       aCamera->SetViewUp(-aViewUp[2], aViewUp[0], aViewUp[1]);
85     
86     // update the camera
87     aCamera->SetPosition(aCenter[0]+aDistance*aViewPlaneNormal[0],
88                          aCenter[1]+aDistance*aViewPlaneNormal[1],
89                          aCenter[2]+aDistance*aViewPlaneNormal[2]);
90
91     // find size of the window
92     int* aWinSize = theRenderer->GetSize();
93     if(aWinSize[0] < aWinSize[1]) 
94       aWidth *= vtkFloatingPointType(aWinSize[1])/vtkFloatingPointType(aWinSize[0]);
95     
96     if(theUsingZeroFocalPoint) 
97       aWidth *= sqrt(2.0);
98     
99     aCamera->SetParallelScale(aWidth/2.0);
100   }
101
102   ResetCameraClippingRange(theRenderer);
103 }
104
105 /*! Compute the bounds of the visible props*/
106 int
107 ComputeVisiblePropBounds(vtkRenderer* theRenderer, 
108                          vtkFloatingPointType theBounds[6])
109 {
110   int aCount = 0;
111   
112   theBounds[0] = theBounds[2] = theBounds[4] = VTK_LARGE_FLOAT;
113   theBounds[1] = theBounds[3] = theBounds[5] = -VTK_LARGE_FLOAT;
114   
115   // loop through all props
116   vtkActorCollection* aCollection = theRenderer->GetActors();
117   aCollection->InitTraversal();
118   while (vtkActor* aProp = aCollection->GetNextActor()) {
119     // if it's invisible, or has no geometry, we can skip the rest 
120     if(aProp->GetVisibility() && aProp->GetMapper()){
121       if(VTKViewer_Actor* anActor = VTKViewer_Actor::SafeDownCast(aProp))
122         if(anActor->IsInfinitive())
123           continue;
124         
125       vtkFloatingPointType *aBounds = aProp->GetBounds();
126       static vtkFloatingPointType MAX_DISTANCE = 0.9*VTK_LARGE_FLOAT;
127       // make sure we haven't got bogus bounds
128       if ( aBounds != NULL &&
129            aBounds[0] > -MAX_DISTANCE && aBounds[1] < MAX_DISTANCE &&
130            aBounds[2] > -MAX_DISTANCE && aBounds[3] < MAX_DISTANCE &&
131            aBounds[4] > -MAX_DISTANCE && aBounds[5] < MAX_DISTANCE )
132       {
133         aCount++;
134
135         theBounds[0] = min(aBounds[0],theBounds[0]);
136         theBounds[2] = min(aBounds[2],theBounds[2]);
137         theBounds[4] = min(aBounds[4],theBounds[4]);
138
139         theBounds[1] = max(aBounds[1],theBounds[1]);
140         theBounds[3] = max(aBounds[3],theBounds[3]);
141         theBounds[5] = max(aBounds[5],theBounds[5]);
142
143       }//not bogus
144     }
145   }
146   return aCount;
147 }
148
149 /*!@see vtkRenderer::ResetCameraClippingRange(vtkFloatingPointType bounds[6]) method*/
150 void
151 ResetCameraClippingRange(vtkRenderer* theRenderer)
152 {
153   if(!theRenderer || !theRenderer->VisibleActorCount()) return;
154   
155   vtkCamera* anActiveCamera = theRenderer->GetActiveCamera();
156   if( anActiveCamera == NULL ){
157     return;
158   }
159   
160   // Find the plane equation for the camera view plane
161   vtkFloatingPointType vn[3];
162   anActiveCamera->GetViewPlaneNormal(vn);
163   vtkFloatingPointType  position[3];
164   anActiveCamera->GetPosition(position);
165   
166   vtkFloatingPointType bounds[6];
167   theRenderer->ComputeVisiblePropBounds(bounds);
168   
169   vtkFloatingPointType center[3];
170   center[0] = (bounds[0] + bounds[1])/2.0;
171   center[1] = (bounds[2] + bounds[3])/2.0;
172   center[2] = (bounds[4] + bounds[5])/2.0;
173   
174   vtkFloatingPointType width = sqrt((bounds[1]-bounds[0])*(bounds[1]-bounds[0]) +
175     (bounds[3]-bounds[2])*(bounds[3]-bounds[2]) +
176     (bounds[5]-bounds[4])*(bounds[5]-bounds[4]));
177   
178   vtkFloatingPointType distance = sqrt((position[0]-center[0])*(position[0]-center[0]) +
179        (position[1]-center[1])*(position[1]-center[1]) +
180        (position[2]-center[2])*(position[2]-center[2]));
181   
182   vtkFloatingPointType range[2] = {distance - width/2.0, distance + width/2.0};
183   
184   // Do not let the range behind the camera throw off the calculation.
185   if (range[0] < 0.0) range[0] = 0.0;
186   
187   anActiveCamera->SetClippingRange( range );
188 }
189
190 /*!Compute trihedron size.*/
191 bool
192 ComputeTrihedronSize( vtkRenderer* theRenderer,
193                       vtkFloatingPointType& theNewSize,
194                       const vtkFloatingPointType theSize, 
195                       const vtkFloatingPointType theSizeInPercents )
196 {
197   // calculating diagonal of visible props of the renderer
198   vtkFloatingPointType bnd[ 6 ];
199   if ( ComputeVisiblePropBounds( theRenderer, bnd ) == 0 )
200   {
201     bnd[ 1 ] = bnd[ 3 ] = bnd[ 5 ] = 100;
202     bnd[ 0 ] = bnd[ 2 ] = bnd[ 4 ] = 0;
203   }
204   vtkFloatingPointType aLength = 0;
205
206   aLength = bnd[ 1 ]-bnd[ 0 ];
207   aLength = max( ( bnd[ 3 ] - bnd[ 2 ] ),aLength );
208   aLength = max( ( bnd[ 5 ] - bnd[ 4 ] ),aLength );
209
210   static vtkFloatingPointType EPS_SIZE = 5.0E-3;
211   theNewSize = aLength * theSizeInPercents / 100.0;
212
213   // if the new trihedron size have sufficient difference, then apply the value
214   return fabs( theNewSize - theSize) > theSize * EPS_SIZE ||
215          fabs( theNewSize-theSize ) > theNewSize * EPS_SIZE;
216 }