VTK  9.2.6
vtkLinearTransform.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Visualization Toolkit
4  Module: vtkLinearTransform.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 =========================================================================*/
36 #ifndef vtkLinearTransform_h
37 #define vtkLinearTransform_h
38 
39 #include "vtkCommonTransformsModule.h" // For export macro
41 
42 class VTKCOMMONTRANSFORMS_EXPORT vtkLinearTransform : public vtkHomogeneousTransform
43 {
44 public:
46  void PrintSelf(ostream& os, vtkIndent indent) override;
47 
52  void TransformNormal(const float in[3], float out[3])
53  {
54  this->Update();
55  this->InternalTransformNormal(in, out);
56  }
57 
62  void TransformNormal(const double in[3], double out[3])
63  {
64  this->Update();
65  this->InternalTransformNormal(in, out);
66  }
67 
72  double* TransformNormal(double x, double y, double z) VTK_SIZEHINT(3)
73  {
74  return this->TransformDoubleNormal(x, y, z);
75  }
76  double* TransformNormal(const double normal[3]) VTK_SIZEHINT(3)
77  {
78  return this->TransformDoubleNormal(normal[0], normal[1], normal[2]);
79  }
80 
82 
86  float* TransformFloatNormal(float x, float y, float z) VTK_SIZEHINT(3)
87  {
88  this->InternalFloatPoint[0] = x;
89  this->InternalFloatPoint[1] = y;
90  this->InternalFloatPoint[2] = z;
91  this->TransformNormal(this->InternalFloatPoint, this->InternalFloatPoint);
92  return this->InternalFloatPoint;
93  }
94  float* TransformFloatNormal(const float normal[3]) VTK_SIZEHINT(3)
95  {
96  return this->TransformFloatNormal(normal[0], normal[1], normal[2]);
97  }
99 
101 
105  double* TransformDoubleNormal(double x, double y, double z) VTK_SIZEHINT(3)
106  {
107  this->InternalDoublePoint[0] = x;
108  this->InternalDoublePoint[1] = y;
109  this->InternalDoublePoint[2] = z;
110  this->TransformNormal(this->InternalDoublePoint, this->InternalDoublePoint);
111  return this->InternalDoublePoint;
112  }
113  double* TransformDoubleNormal(const double normal[3]) VTK_SIZEHINT(3)
114  {
115  return this->TransformDoubleNormal(normal[0], normal[1], normal[2]);
116  }
118 
123  double* TransformVector(double x, double y, double z) VTK_SIZEHINT(3)
124  {
125  return this->TransformDoubleVector(x, y, z);
126  }
127  double* TransformVector(const double normal[3]) VTK_SIZEHINT(3)
128  {
129  return this->TransformDoubleVector(normal[0], normal[1], normal[2]);
130  }
131 
136  void TransformVector(const float in[3], float out[3])
137  {
138  this->Update();
139  this->InternalTransformVector(in, out);
140  }
141 
146  void TransformVector(const double in[3], double out[3])
147  {
148  this->Update();
149  this->InternalTransformVector(in, out);
150  }
151 
153 
157  float* TransformFloatVector(float x, float y, float z) VTK_SIZEHINT(3)
158  {
159  this->InternalFloatPoint[0] = x;
160  this->InternalFloatPoint[1] = y;
161  this->InternalFloatPoint[2] = z;
162  this->TransformVector(this->InternalFloatPoint, this->InternalFloatPoint);
163  return this->InternalFloatPoint;
164  }
165  float* TransformFloatVector(const float vec[3]) VTK_SIZEHINT(3)
166  {
167  return this->TransformFloatVector(vec[0], vec[1], vec[2]);
168  }
170 
172 
176  double* TransformDoubleVector(double x, double y, double z) VTK_SIZEHINT(3)
177  {
178  this->InternalDoublePoint[0] = x;
179  this->InternalDoublePoint[1] = y;
180  this->InternalDoublePoint[2] = z;
181  this->TransformVector(this->InternalDoublePoint, this->InternalDoublePoint);
182  return this->InternalDoublePoint;
183  }
184  double* TransformDoubleVector(const double vec[3]) VTK_SIZEHINT(3)
185  {
186  return this->TransformDoubleVector(vec[0], vec[1], vec[2]);
187  }
189 
194  void TransformPoints(vtkPoints* inPts, vtkPoints* outPts) override;
195 
200  virtual void TransformNormals(vtkDataArray* inNms, vtkDataArray* outNms);
201 
206  virtual void TransformVectors(vtkDataArray* inVrs, vtkDataArray* outVrs);
207 
212  void TransformPointsNormalsVectors(vtkPoints* inPts, vtkPoints* outPts, vtkDataArray* inNms,
213  vtkDataArray* outNms, vtkDataArray* inVrs, vtkDataArray* outVrs, int nOptionalVectors = 0,
214  vtkDataArray** inVrsArr = nullptr, vtkDataArray** outVrsArr = nullptr) override;
215 
221  {
222  return static_cast<vtkLinearTransform*>(this->GetInverse());
223  }
224 
226 
230  void InternalTransformPoint(const float in[3], float out[3]) override;
231  void InternalTransformPoint(const double in[3], double out[3]) override;
233 
235 
239  virtual void InternalTransformNormal(const float in[3], float out[3]);
240  virtual void InternalTransformNormal(const double in[3], double out[3]);
242 
244 
248  virtual void InternalTransformVector(const float in[3], float out[3]);
249  virtual void InternalTransformVector(const double in[3], double out[3]);
251 
253 
259  const float in[3], float out[3], float derivative[3][3]) override;
261  const double in[3], double out[3], double derivative[3][3]) override;
263 
264 protected:
265  vtkLinearTransform() = default;
266  ~vtkLinearTransform() override = default;
267 
268 private:
269  vtkLinearTransform(const vtkLinearTransform&) = delete;
270  void operator=(const vtkLinearTransform&) = delete;
271 };
272 
273 #endif
void TransformNormal(const double in[3], double out[3])
Apply the transformation to a double-precision normal.
double * TransformDoubleNormal(double x, double y, double z)
Apply the transformation to a double-precision (x,y,z) normal.
void TransformVector(const float in[3], float out[3])
Apply the transformation to a vector.
float * TransformFloatNormal(float x, float y, float z)
Apply the transformation to an (x,y,z) normal.
superclass for homogeneous transformations
void InternalTransformPoint(const float in[3], float out[3]) override
This will calculate the transformation without calling Update.
double * TransformDoubleNormal(const double normal[3])
Apply the transformation to a double-precision (x,y,z) normal.
float * TransformFloatNormal(const float normal[3])
Apply the transformation to an (x,y,z) normal.
vtkAbstractTransform * GetInverse()
Get the inverse of this transform.
double * TransformVector(double x, double y, double z)
Synonymous with TransformDoubleVector(x,y,z).
void TransformVector(const double in[3], double out[3])
Apply the transformation to a double-precision vector.
double * TransformVector(const double normal[3])
vtkLinearTransform * GetLinearInverse()
Just like GetInverse, but it includes a typecast to vtkLinearTransform.
float * TransformFloatVector(const float vec[3])
Apply the transformation to an (x,y,z) vector.
void Update()
Update the transform to account for any changes which have been made.
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
a simple class to control print indentation
Definition: vtkIndent.h:39
double * TransformNormal(const double normal[3])
abstract superclass for arrays of numeric data
Definition: vtkDataArray.h:55
#define VTK_SIZEHINT(...)
double * TransformNormal(double x, double y, double z)
Synonymous with TransformDoubleNormal(x,y,z).
void TransformNormal(const float in[3], float out[3])
Apply the transformation to a normal.
void TransformPoints(vtkPoints *inPts, vtkPoints *outPts) override
Apply the transformation to a series of points, and append the results to outPts. ...
void InternalTransformDerivative(const float in[3], float out[3], float derivative[3][3]) override
This will calculate the transformation as well as its derivative without calling Update.
float * TransformFloatVector(float x, float y, float z)
Apply the transformation to an (x,y,z) vector.
double * TransformDoubleVector(const double vec[3])
Apply the transformation to a double-precision (x,y,z) vector.
void TransformPointsNormalsVectors(vtkPoints *inPts, vtkPoints *outPts, vtkDataArray *inNms, vtkDataArray *outNms, vtkDataArray *inVrs, vtkDataArray *outVrs, int nOptionalVectors=0, vtkDataArray **inVrsArr=nullptr, vtkDataArray **outVrsArr=nullptr) override
Apply the transformation to a combination of points, normals and vectors.
double * TransformDoubleVector(double x, double y, double z)
Apply the transformation to a double-precision (x,y,z) vector.
represent and manipulate 3D points
Definition: vtkPoints.h:39
abstract superclass for linear transformations