VTK  9.3.20240418
vtkTransform2D.h
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2 // SPDX-License-Identifier: BSD-3-Clause
3 
29 #ifndef vtkTransform2D_h
30 #define vtkTransform2D_h
31 
32 #include "vtkCommonTransformsModule.h" // For export macro
33 #include "vtkObject.h"
34 
35 #include "vtkMatrix3x3.h" // Needed for inline methods
36 
37 VTK_ABI_NAMESPACE_BEGIN
38 class vtkPoints2D;
39 
40 class VTKCOMMONTRANSFORMS_EXPORT vtkTransform2D : public vtkObject
41 {
42 public:
43  static vtkTransform2D* New();
44  vtkTypeMacro(vtkTransform2D, vtkObject);
45  void PrintSelf(ostream& os, vtkIndent indent) override;
46 
50  void Identity();
51 
55  void Inverse();
56 
61  void Translate(double x, double y);
62  void Translate(const double x[2]) { this->Translate(x[0], x[1]); }
63  void Translate(const float x[2]) { this->Translate(x[0], x[1]); }
64 
69  void Rotate(double angle);
70 
75  void Scale(double x, double y);
76  void Scale(const double s[2]) { this->Scale(s[0], s[1]); }
77  void Scale(const float s[2]) { this->Scale(s[0], s[1]); }
78 
82  void SetMatrix(vtkMatrix3x3* matrix) { this->SetMatrix(matrix->GetData()); }
83  void SetMatrix(const double elements[9]);
84 
86 
89  vtkGetObjectMacro(Matrix, vtkMatrix3x3);
90  void GetMatrix(vtkMatrix3x3* matrix);
92 
94 
99  void GetPosition(double pos[2]);
100  void GetPosition(float pos[2])
101  {
102  double temp[2];
103  this->GetPosition(temp);
104  pos[0] = static_cast<float>(temp[0]);
105  pos[1] = static_cast<float>(temp[1]);
106  }
108 
110 
115  void GetScale(double scale[2]);
116  void GetScale(float pos[2])
117  {
118  double temp[2];
119  this->GetScale(temp);
120  pos[0] = static_cast<float>(temp[0]);
121  pos[1] = static_cast<float>(temp[1]);
122  }
124 
129  void GetInverse(vtkMatrix3x3* inverse);
130 
136  void GetTranspose(vtkMatrix3x3* transpose);
137 
141  vtkMTimeType GetMTime() override;
142 
148  void TransformPoints(const float* inPts, float* outPts, int n);
149 
155  void TransformPoints(const double* inPts, double* outPts, int n);
156 
161  void TransformPoints(vtkPoints2D* inPts, vtkPoints2D* outPts);
162 
168  void InverseTransformPoints(const float* inPts, float* outPts, int n);
169 
175  void InverseTransformPoints(const double* inPts, double* outPts, int n);
176 
182 
184 
189  void MultiplyPoint(const float in[3], float out[3]) { this->GetMatrix()->MultiplyPoint(in, out); }
190  void MultiplyPoint(const double in[3], double out[3])
191  {
192  this->GetMatrix()->MultiplyPoint(in, out);
193  }
195 
196 protected:
198  ~vtkTransform2D() override;
199 
201 
204 
205 private:
206  vtkTransform2D(const vtkTransform2D&) = delete;
207  void operator=(const vtkTransform2D&) = delete;
208 };
209 
210 VTK_ABI_NAMESPACE_END
211 #endif
a simple class to control print indentation
Definition: vtkIndent.h:108
represent and manipulate 3x3 transformation matrices
Definition: vtkMatrix3x3.h:56
double * GetData()
Return a pointer to the first element of the matrix (double[9]).
Definition: vtkMatrix3x3.h:199
abstract base class for most VTK objects
Definition: vtkObject.h:162
represent and manipulate 2D points
Definition: vtkPoints2D.h:27
describes linear transformations via a 3x3 matrix
vtkMTimeType GetMTime() override
Override GetMTime to account for input and concatenation.
void SetMatrix(vtkMatrix3x3 *matrix)
Set the current matrix directly.
void GetScale(float pos[2])
Return the x and y scale from the current transformation matrix as an array of two floating point num...
void Translate(const float x[2])
void MultiplyPoint(const double in[3], double out[3])
Use this method only if you wish to compute the transformation in homogeneous (x,y,...
void GetTranspose(vtkMatrix3x3 *transpose)
Return a matrix which is the transpose of the current transformation matrix.
void GetPosition(float pos[2])
Return the position from the current transformation matrix as an array of two floating point numbers.
void Rotate(double angle)
Create a rotation matrix and concatenate it with the current transformation.
vtkMatrix3x3 * InverseMatrix
void Scale(const float s[2])
void TransformPoints(const float *inPts, float *outPts, int n)
Apply the transformation to a series of points, and append the results to outPts.
void SetMatrix(const double elements[9])
vtkMatrix3x3 * Matrix
void GetInverse(vtkMatrix3x3 *inverse)
Return a matrix which is the inverse of the current transformation matrix.
void InverseTransformPoints(const float *inPts, float *outPts, int n)
Apply the transformation to a series of points, and append the results to outPts.
~vtkTransform2D() override
void MultiplyPoint(const float in[3], float out[3])
Use this method only if you wish to compute the transformation in homogeneous (x,y,...
void Scale(double x, double y)
Create a scale matrix (i.e.
void Inverse()
Invert the transformation.
void InternalDeepCopy(vtkTransform2D *t)
static vtkTransform2D * New()
void PrintSelf(ostream &os, vtkIndent indent) override
Methods invoked by print to print information about the object including superclasses.
void Identity()
Set the transformation to the identity transformation.
void GetScale(double scale[2])
Return the x and y scale from the current transformation matrix as an array of two floating point num...
void Scale(const double s[2])
void GetMatrix(vtkMatrix3x3 *matrix)
Get the underlying 3x3 matrix.
void InverseTransformPoints(vtkPoints2D *inPts, vtkPoints2D *outPts)
Apply the transformation to a series of points, and append the results to outPts.
void InverseTransformPoints(const double *inPts, double *outPts, int n)
Apply the transformation to a series of points, and append the results to outPts.
void TransformPoints(const double *inPts, double *outPts, int n)
Apply the transformation to a series of points, and append the results to outPts.
void Translate(const double x[2])
void GetPosition(double pos[2])
Return the position from the current transformation matrix as an array of two floating point numbers.
void Translate(double x, double y)
Create a translation matrix and concatenate it with the current transformation.
void TransformPoints(vtkPoints2D *inPts, vtkPoints2D *outPts)
Apply the transformation to a series of points, and append the results to outPts.
@ scale
Definition: vtkX3D.h:229
vtkTypeUInt32 vtkMTimeType
Definition: vtkType.h:270