ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkRigid3DTransform.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright Insight Software Consortium
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkRigid3DTransform_h
19 #define itkRigid3DTransform_h
20 
21 #include <iostream>
23 #include "itkVersor.h"
24 
25 namespace itk
26 {
56 template<typename TParametersValueType=double>
58  public MatrixOffsetTransformBase<TParametersValueType, 3, 3>
59 {
60 public:
66 
67 #ifdef ITKV3_COMPATIBILITY
68 
69  itkNewMacro(Self);
70 #endif
71 
74 
76  itkStaticConstMacro(SpaceDimension, unsigned int, 3);
77  itkStaticConstMacro(InputSpaceDimension, unsigned int, 3);
78  itkStaticConstMacro(OutputSpaceDimension, unsigned int, 3);
79  itkStaticConstMacro(ParametersDimension, unsigned int, 12);
81 
103 
107  typedef typename InverseTransformBaseType::Pointer InverseTransformBasePointer;
108 
119  virtual void SetParameters(const ParametersType & parameters) ITK_OVERRIDE;
120 
126  virtual void SetMatrix(const MatrixType & matrix) ITK_OVERRIDE;
127 
133  virtual void SetMatrix(const MatrixType & matrix, const TParametersValueType tolerance );
134 
142  void Translate(const OffsetType & offset, bool pre = false);
143 
148  bool MatrixIsOrthogonal(const MatrixType & matrix,
149  const TParametersValueType tolerance =
151 
152 #ifdef ITKV3_COMPATIBILITY
153 
154  //NOTE: itkLegacyRemove can not be used for GetInverse
155  // because in itkV3 mode these functions
156  // must be traversed when calling the child classes
157  // member functions
158  // (with no real effect) for backwards compatibility.
159  // In ITKv4 mode only the super class is needed
160  bool GetInverse(Self *inverse) const;
161 
163  //NOTE: itkLegacyRemove can not be used for GetInverseTransform
164  // because in itkV3 mode these functions
165  // must be traversed when calling the child classes
166  // member functions
167  // (with no real effect) for backwards compatibility.
168  // In ITKv4 mode only the super class is needed
170 
179  itkLegacyMacro(const MatrixType & GetRotationMatrix() const);
180 
191  itkLegacyMacro(virtual void SetRotationMatrix(const MatrixType & matrix) );
192 #endif
193 
205  itkLegacyMacro(InputPointType BackTransform(const OutputPointType & point) const);
206  itkLegacyMacro(InputVectorType BackTransform(const OutputVectorType & vector) const);
207  itkLegacyMacro(InputVnlVectorType BackTransform(const OutputVnlVectorType & vector) const);
208  itkLegacyMacro(InputCovariantVectorType BackTransform(const OutputCovariantVectorType & vector) const);
210 
211 protected:
212  Rigid3DTransform(const MatrixType & matrix,
213  const OutputVectorType & offset);
214  Rigid3DTransform(unsigned int paramDim);
217 
221  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
222 
223 private:
224  Rigid3DTransform(const Self &) ITK_DELETE_FUNCTION;
225  void operator=(const Self &) ITK_DELETE_FUNCTION;
226 }; //class Rigid3DTransform
227 } // namespace itk
228 
229 #ifndef ITK_MANUAL_INSTANTIATION
230 #include "itkRigid3DTransform.hxx"
231 #endif
232 
233 #endif /* itkRigid3DTransform_h */
Superclass::InputCovariantVectorType InputCovariantVectorType
Superclass::FixedParametersType FixedParametersType
static const unsigned int OutputSpaceDimension
Point< TParametersValueType, itkGetStaticConstMacro(OutputSpaceDimension)> OutputPointType
CovariantVector< TParametersValueType, itkGetStaticConstMacro(InputSpaceDimension)> InputCovariantVectorType
Matrix< TParametersValueType, itkGetStaticConstMacro(OutputSpaceDimension), itkGetStaticConstMacro(InputSpaceDimension)> MatrixType
Superclass::TranslationType TranslationType
Superclass::InverseMatrixType InverseMatrixType
Light weight base class for most itk classes.
Superclass::CenterType CenterType
Rigid3DTransform of a vector space (e.g. space coordinates)
CovariantVector< TParametersValueType, itkGetStaticConstMacro(OutputSpaceDimension)> OutputCovariantVectorType
Matrix and Offset transformation of a vector space (e.g. space coordinates)
virtual void SetParameters(const ParametersType &parameters) override
InverseTransformBaseType::Pointer InverseTransformBasePointer
Superclass::MatrixType MatrixType
Superclass::OutputVectorType OutputVectorType
itkLegacyMacro(InputPointType BackTransform(const OutputPointType &point) const)
vnl_vector_fixed< TParametersValueType, itkGetStaticConstMacro(OutputSpaceDimension)> OutputVnlVectorType
static const unsigned int ParametersDimension
Superclass::InputVectorType InputVectorType
void Translate(const OffsetType &offset, bool pre=false)
Superclass::OutputVectorValueType OutputVectorValueType
Matrix< TParametersValueType, itkGetStaticConstMacro(InputSpaceDimension), itkGetStaticConstMacro(OutputSpaceDimension)> InverseMatrixType
bool MatrixIsOrthogonal(const MatrixType &matrix, const TParametersValueType tolerance=MatrixOrthogonalityTolerance< TParametersValueType >::GetTolerance())
Superclass::OutputPointType OutputPointType
Superclass::OutputCovariantVectorType OutputCovariantVectorType
Superclass::JacobianType JacobianType
SmartPointer< const Self > ConstPointer
Superclass::ParametersValueType ParametersValueType
SmartPointer< Self > Pointer
virtual void SetMatrix(const MatrixType &matrix) override
Superclass::InputVnlVectorType InputVnlVectorType
void PrintSelf(std::ostream &os, Indent indent) const override
Vector< TParametersValueType, itkGetStaticConstMacro(OutputSpaceDimension)> OutputVectorType
static const unsigned int InputSpaceDimension
Superclass::OffsetType OffsetType
Superclass::ParametersType ParametersType
MatrixOffsetTransformBase< TParametersValueType, 3, 3 > Superclass
Control indentation during Print() invocation.
Definition: itkIndent.h:49
vnl_vector_fixed< TParametersValueType, itkGetStaticConstMacro(InputSpaceDimension)> InputVnlVectorType
static const unsigned int SpaceDimension
virtual InverseTransformBasePointer GetInverseTransform() const override
Superclass::ScalarType ScalarType
Superclass::InverseTransformBaseType InverseTransformBaseType
Point< TParametersValueType, itkGetStaticConstMacro(InputSpaceDimension)> InputPointType
Superclass::InputPointType InputPointType
Superclass::OutputVnlVectorType OutputVnlVectorType
Superclass::MatrixValueType MatrixValueType
Vector< TParametersValueType, itkGetStaticConstMacro(InputSpaceDimension)> InputVectorType
Superclass::FixedParametersValueType FixedParametersValueType