ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkSimilarity3DTransform.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 itkSimilarity3DTransform_h
19 #define itkSimilarity3DTransform_h
20 
21 #include <iostream>
23 
24 namespace itk
25 {
46 template<typename TParametersValueType=double>
48  public VersorRigid3DTransform<TParametersValueType>
49 {
50 public:
56 
58  itkNewMacro(Self);
59 
62 
64  itkStaticConstMacro(SpaceDimension, unsigned int, 3);
65  itkStaticConstMacro(InputSpaceDimension, unsigned int, 3);
66  itkStaticConstMacro(OutputSpaceDimension, unsigned int, 3);
67  itkStaticConstMacro(ParametersDimension, unsigned int, 7);
69 
88 
91  typedef typename Superclass::AxisType AxisType;
92  typedef typename Superclass::AngleType AngleType;
93  typedef TParametersValueType ScaleType;
94 
96  virtual void SetIdentity(void) ITK_OVERRIDE;
97 
104  virtual void SetMatrix(const MatrixType & matrix) ITK_OVERRIDE;
105 
112  virtual void SetMatrix(const MatrixType & matrix, double tolerance) ITK_OVERRIDE;
113 
118  void SetParameters(const ParametersType & parameters) ITK_OVERRIDE;
119 
120  virtual const ParametersType & GetParameters(void) const ITK_OVERRIDE;
121 
123  void SetScale(ScaleType scale);
124 
125  itkGetConstReferenceMacro(Scale, ScaleType);
126 
131  virtual void ComputeJacobianWithRespectToParameters( const InputPointType & p, JacobianType & jacobian) const ITK_OVERRIDE;
132 
133 protected:
134  Similarity3DTransform(const MatrixType & matrix, const OutputVectorType & offset);
135  Similarity3DTransform(unsigned int paramDim);
138  {
139  }
140 
141  void PrintSelf(std::ostream & os, Indent indent) const ITK_OVERRIDE;
142 
145  void ComputeMatrix() ITK_OVERRIDE;
146 
148  void ComputeMatrixParameters() ITK_OVERRIDE;
149 
150 private:
151  Similarity3DTransform(const Self &); // purposely not implemented
152  void operator=(const Self &); // purposely not implemented
153 
155 }; // class Similarity3DTransform
156 } // namespace itk
157 
158 #ifndef ITK_MANUAL_INSTANTIATION
159 #include "itkSimilarity3DTransform.hxx"
160 #endif
161 
162 #endif /* itkSimilarity3DTransform_h */
void ComputeMatrixParameters() override
Superclass::InputCovariantVectorType InputCovariantVectorType
Superclass::OutputVnlVectorType OutputVnlVectorType
Superclass::InverseMatrixType InverseMatrixType
Superclass::InputPointType InputPointType
Superclass::OutputPointType OutputPointType
Light weight base class for most itk classes.
Superclass::InputVectorType InputVectorType
Superclass::FixedParametersType FixedParametersType
Superclass::InputVnlVectorType InputVnlVectorType
Superclass::MatrixType MatrixType
VersorRigid3DTransform of a vector space (e.g. space coordinates)
Superclass::InputCovariantVectorType InputCovariantVectorType
Superclass::InputVnlVectorType InputVnlVectorType
static const unsigned int OutputSpaceDimension
void SetParameters(const ParametersType &parameters) override
static const unsigned int ParametersDimension
VersorRigid3DTransform< TParametersValueType > Superclass
Superclass::MatrixType MatrixType
Superclass::OutputVectorType OutputVectorType
Superclass::InputVectorType InputVectorType
Superclass::TranslationType TranslationType
void ComputeMatrix() override
void PrintSelf(std::ostream &os, Indent indent) const override
virtual const ParametersType & GetParameters(void) const override
Superclass::FixedParametersType FixedParametersType
Superclass::OutputVectorType OutputVectorType
Superclass::OutputVectorType OutputVectorType
Superclass::ScalarType ScalarType
void SetScale(ScaleType scale)
Similarity3DTransform of a vector space (e.g. space coordinates)
Superclass::JacobianType JacobianType
Superclass::VersorType VersorType
Superclass::OutputPointType OutputPointType
Superclass::ParametersType ParametersType
virtual void SetMatrix(const MatrixType &matrix) override
Superclass::JacobianType JacobianType
Superclass::OffsetType OffsetType
Superclass::OutputVnlVectorType OutputVnlVectorType
Superclass::OutputCovariantVectorType OutputCovariantVectorType
Superclass::InputPointType InputPointType
Superclass::JacobianType JacobianType
Superclass::TranslationType TranslationType
virtual void ComputeJacobianWithRespectToParameters(const InputPointType &p, JacobianType &jacobian) const override
Superclass::InverseMatrixType InverseMatrixType
Superclass::ParametersType ParametersType
Superclass::ParametersType ParametersType
static const unsigned int InputSpaceDimension
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Superclass::CenterType CenterType
virtual void SetIdentity(void) override
Superclass::InputPointType InputPointType
Superclass::OutputCovariantVectorType OutputCovariantVectorType
SmartPointer< const Self > ConstPointer
static const unsigned int SpaceDimension