VTK  9.0.1
vtkBlockSortHelper.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkBlockSortHelper.h
5 
6  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7  All rights reserved.
8  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9 
10  This software is distributed WITHOUT ANY WARRANTY; without even
11  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12  PURPOSE. See the above copyright notice for more information.
13 
14 =========================================================================*/
20 #ifndef vtkBlockSortHelper_h
21 #define vtkBlockSortHelper_h
22 
23 #include "vtkCamera.h"
24 #include "vtkImageData.h"
25 #include "vtkMatrix4x4.h"
26 #include "vtkNew.h"
27 #include "vtkRenderer.h"
28 #include "vtkVolumeMapper.h"
29 
31 {
32 
39 template <typename T>
41 {
42  double CameraPosition[4];
43 
44  //----------------------------------------------------------------------------
46  {
47  vtkCamera* cam = ren->GetActiveCamera();
48  double camWorldPos[4];
49 
50  cam->GetPosition(camWorldPos);
51  camWorldPos[3] = 1.0;
52 
53  // Transform the camera position to the volume (dataset) coordinate system.
54  vtkNew<vtkMatrix4x4> InverseVolumeMatrix;
55  InverseVolumeMatrix->DeepCopy(volMatrix);
56  InverseVolumeMatrix->Invert();
57  InverseVolumeMatrix->MultiplyPoint(camWorldPos, this->CameraPosition);
58  };
59 
60  //----------------------------------------------------------------------------
61  bool operator()(T* first, T* second);
62 
72  //----------------------------------------------------------------------------
74  {
75  double center[3];
76  double bounds[6];
77 
78  first->GetBounds(bounds);
79  this->ComputeCenter(bounds, center);
80  double const dist1 = vtkMath::Distance2BetweenPoints(center, this->CameraPosition);
81 
82  second->GetBounds(bounds);
83  this->ComputeCenter(bounds, center);
84  double const dist2 = vtkMath::Distance2BetweenPoints(center, this->CameraPosition);
85 
86  return dist2 < dist1;
87  };
88 
89  //----------------------------------------------------------------------------
90  inline void ComputeCenter(double const* bounds, double* center)
91  {
92  center[0] = bounds[0] + std::abs(bounds[1] - bounds[0]) / 2.0;
93  center[1] = bounds[2] + std::abs(bounds[3] - bounds[2]) / 2.0;
94  center[2] = bounds[4] + std::abs(bounds[5] - bounds[4]) / 2.0;
95  };
96 };
97 
98 //----------------------------------------------------------------------------
99 template <>
101 {
102  return CompareByDistanceDescending(first, second);
103 };
104 
105 //----------------------------------------------------------------------------
106 template <>
108  vtkVolumeMapper* first, vtkVolumeMapper* second)
109 {
110  vtkImageData* firstIm = first->GetInput();
111  vtkImageData* secondIm = second->GetInput();
112 
113  return CompareByDistanceDescending(firstIm, secondIm);
114 };
115 }
116 
117 #endif // vtkBlockSortHelper_h
118 // VTK-HeaderTest-Exclude: vtkBlockSortHelper.h
bool operator()(T *first, T *second)
Abstract class for a volume mapper.
represent and manipulate 4x4 transformation matrices
Definition: vtkMatrix4x4.h:35
double * GetBounds()
Return a pointer to the geometry bounding box in the form (xmin,xmax, ymin,ymax, zmin,zmax).
abstract specification for renderers
Definition: vtkRenderer.h:67
vtkCamera * GetActiveCamera()
Get the current camera.
virtual double * GetPosition()
Set/Get the position of the camera in world coordinates.
void MultiplyPoint(const float in[4], float out[4])
Multiply a homogeneous coordinate by this matrix, i.e.
Definition: vtkMatrix4x4.h:136
void DeepCopy(const vtkMatrix4x4 *source)
Set the elements of the matrix to the same values as the elements of the given source matrix...
Definition: vtkMatrix4x4.h:53
BackToFront(vtkRenderer *ren, vtkMatrix4x4 *volMatrix)
static void Invert(const vtkMatrix4x4 *in, vtkMatrix4x4 *out)
Matrix Inversion (adapted from Richard Carling in "Graphics Gems," Academic Press, 1990).
Definition: vtkMatrix4x4.h:113
a virtual camera for 3D rendering
Definition: vtkCamera.h:45
operator() for back-to-front sorting.
topologically and geometrically regular array of data
Definition: vtkImageData.h:41
virtual vtkImageData * GetInput()
Set/Get the input data.
Collection of comparison functions for std::sort.
static float Distance2BetweenPoints(const float p1[3], const float p2[3])
Compute distance squared between two points p1 and p2.
Definition: vtkMath.h:1426
void ComputeCenter(double const *bounds, double *center)
bool CompareByDistanceDescending(vtkImageData *first, vtkImageData *second)
Compares distances from images (first, second) to the camera position.