ITK  5.2.0
Insight Toolkit
itkMatrixOffsetTransformBase.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
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 itkMatrixOffsetTransformBase_h
19 #define itkMatrixOffsetTransformBase_h
20 
21 
22 #include "itkMacro.h"
23 #include "itkMatrix.h"
24 #include "itkTransform.h"
25 
26 #include <iostream>
27 
28 namespace itk
29 {
30 
31 /* MatrixOrthogonalityTolerance is a utility to
32  * allow setting the tolerance limits used for
33  * checking if a matrix meet the orthogonality
34  * constraints of being a rigid rotation matrix.
35  * The tolerance needs to be different for
36  * matrices of type float vs. double.
37  */
38 template <typename T>
40 
41 template <>
42 class ITK_TEMPLATE_EXPORT MatrixOrthogonalityTolerance<double>
43 {
44 public:
45  static double
47  {
48  return 1e-10;
49  }
50 };
51 
52 template <>
53 class ITK_TEMPLATE_EXPORT MatrixOrthogonalityTolerance<float>
54 {
55 public:
56  static float
58  {
59  return 1e-5f;
60  }
61 };
62 
105 template <typename TParametersValueType = double, unsigned int NInputDimensions = 3, unsigned int NOutputDimensions = 3>
106 class ITK_TEMPLATE_EXPORT MatrixOffsetTransformBase
107  : public Transform<TParametersValueType, NInputDimensions, NOutputDimensions>
108 {
109 public:
110  ITK_DISALLOW_COPY_AND_MOVE(MatrixOffsetTransformBase);
111 
115 
118 
120  itkTypeMacro(MatrixOffsetTransformBase, Transform);
121 
123  itkNewMacro(Self);
124 
126  static constexpr unsigned int InputSpaceDimension = NInputDimensions;
127  static constexpr unsigned int OutputSpaceDimension = NOutputDimensions;
128  static constexpr unsigned int ParametersDimension = NOutputDimensions * (NInputDimensions + 1);
129 
131  using FixedParametersType = typename Superclass::FixedParametersType;
132  using FixedParametersValueType = typename Superclass::FixedParametersValueType;
133  using ParametersType = typename Superclass::ParametersType;
134  using ParametersValueType = typename Superclass::ParametersValueType;
135 
137  using JacobianType = typename Superclass::JacobianType;
138  using JacobianPositionType = typename Superclass::JacobianPositionType;
139  using InverseJacobianPositionType = typename Superclass::InverseJacobianPositionType;
140 
142  using TransformCategoryEnum = typename Superclass::TransformCategoryEnum;
143 
145  using ScalarType = typename Superclass::ScalarType;
146 
151 
155 
156  using InputVectorPixelType = typename Superclass::InputVectorPixelType;
157  using OutputVectorPixelType = typename Superclass::OutputVectorPixelType;
158 
160  using InputDiffusionTensor3DType = typename Superclass::InputDiffusionTensor3DType;
161  using OutputDiffusionTensor3DType = typename Superclass::OutputDiffusionTensor3DType;
162 
164  using InputSymmetricSecondRankTensorType = typename Superclass::InputSymmetricSecondRankTensorType;
165  using OutputSymmetricSecondRankTensorType = typename Superclass::OutputSymmetricSecondRankTensorType;
166 
168 
170  using InputVnlVectorType = vnl_vector_fixed<TParametersValueType, Self::InputSpaceDimension>;
171  using OutputVnlVectorType = vnl_vector_fixed<TParametersValueType, Self::OutputSpaceDimension>;
172 
178 
182 
185 
187 
190 
192 
194 
197  using InverseTransformBaseType = typename Superclass::InverseTransformBaseType;
199 
203  virtual void
204  SetIdentity();
205 
210  GetTransformCategory() const override
211  {
212  return Self::TransformCategoryEnum::Linear;
213  }
214 
226  virtual void
227  SetMatrix(const MatrixType & matrix)
228  {
229  m_Matrix = matrix;
230  this->ComputeOffset();
231  this->ComputeMatrixParameters();
232  m_MatrixMTime.Modified();
233  this->Modified();
234  return;
235  }
237 
245  virtual const MatrixType &
246  GetMatrix() const
247  {
248  return m_Matrix;
249  }
250 
259  void
260  SetOffset(const OutputVectorType & offset)
261  {
262  m_Offset = offset;
263  this->ComputeTranslation();
264  this->Modified();
265  return;
266  }
268 
274  const OutputVectorType &
275  GetOffset() const
276  {
277  return m_Offset;
278  }
279 
302  void
303  SetCenter(const InputPointType & center)
304  {
305  m_Center = center;
306  this->ComputeOffset();
307  this->Modified();
308  return;
309  }
311 
318  const InputPointType &
319  GetCenter() const
320  {
321  return m_Center;
322  }
323 
330  void
331  SetTranslation(const OutputVectorType & translation)
332  {
333  m_Translation = translation;
334  this->ComputeOffset();
335  this->Modified();
336  return;
337  }
339 
346  const OutputVectorType &
348  {
349  return m_Translation;
350  }
351 
356  void
357  SetParameters(const ParametersType & parameters) override;
358 
360  const ParametersType &
361  GetParameters() const override;
362 
364  void
365  SetFixedParameters(const FixedParametersType &) override;
366 
368  const FixedParametersType &
369  GetFixedParameters() const override;
370 
382  void
383  Compose(const Self * other, bool pre = false);
384 
393  OutputPointType
394  TransformPoint(const InputPointType & point) const override;
395 
396  using Superclass::TransformVector;
397 
398  OutputVectorType
399  TransformVector(const InputVectorType & vect) const override;
400 
401  OutputVnlVectorType
402  TransformVector(const InputVnlVectorType & vect) const override;
403 
404  OutputVectorPixelType
405  TransformVector(const InputVectorPixelType & vect) const override;
406 
407  using Superclass::TransformCovariantVector;
408 
409  OutputCovariantVectorType
410  TransformCovariantVector(const InputCovariantVectorType & vec) const override;
411 
412  OutputVectorPixelType
413  TransformCovariantVector(const InputVectorPixelType & vect) const override;
414 
415  using Superclass::TransformDiffusionTensor3D;
416 
417  OutputDiffusionTensor3DType
418  TransformDiffusionTensor3D(const InputDiffusionTensor3DType & tensor) const override;
419 
420  OutputVectorPixelType
421  TransformDiffusionTensor3D(const InputVectorPixelType & tensor) const override;
422 
423  using Superclass::TransformSymmetricSecondRankTensor;
424  OutputSymmetricSecondRankTensorType
425  TransformSymmetricSecondRankTensor(const InputSymmetricSecondRankTensorType & inputTensor) const override;
426 
427  OutputVectorPixelType
428  TransformSymmetricSecondRankTensor(const InputVectorPixelType & inputTensor) const override;
429 
430 
440  void
441  ComputeJacobianWithRespectToParameters(const InputPointType & p, JacobianType & jacobian) const override;
442 
443 
447  void
448  ComputeJacobianWithRespectToPosition(const InputPointType & x, JacobianPositionType & jac) const override;
449  using Superclass::ComputeJacobianWithRespectToPosition;
450 
454  void
455  ComputeInverseJacobianWithRespectToPosition(const InputPointType & x,
456  InverseJacobianPositionType & jac) const override;
457  using Superclass::ComputeInverseJacobianWithRespectToPosition;
458 
477  bool
478  GetInverse(Self * inverse) const;
480 
482  InverseTransformBasePointer
483  GetInverseTransform() const override;
484 
490  bool
491  IsLinear() const override
492  {
493  return true;
494  }
495 
496 protected:
499  const InverseMatrixType &
500  GetInverseMatrix() const;
501 
502 protected:
510  MatrixOffsetTransformBase(const MatrixType & matrix, const OutputVectorType & offset);
511  MatrixOffsetTransformBase(unsigned int paramDims);
514 
516  ~MatrixOffsetTransformBase() override = default;
517 
519  void
520  PrintSelf(std::ostream & os, Indent indent) const override;
521 
522  const InverseMatrixType &
524  {
525  return m_InverseMatrix;
526  }
527  void
529  {
530  m_InverseMatrix = matrix;
531  m_InverseMatrixMTime.Modified();
532  }
533  bool
535  {
536  if (m_MatrixMTime != m_InverseMatrixMTime)
537  {
538  return true;
539  }
540  else
541  {
542  return false;
543  }
544  }
545 
546  virtual void
547  ComputeMatrixParameters();
548 
549  virtual void
550  ComputeMatrix();
551 
552  void
553  SetVarMatrix(const MatrixType & matrix)
554  {
555  m_Matrix = matrix;
556  m_MatrixMTime.Modified();
557  }
558 
559  virtual void
560  ComputeTranslation();
561 
562  void
563  SetVarTranslation(const OutputVectorType & translation)
564  {
565  m_Translation = translation;
566  }
567 
568  virtual void
569  ComputeOffset();
570 
571  void
573  {
574  m_Offset = offset;
575  }
576 
577  void
578  SetVarCenter(const InputPointType & center)
579  {
580  m_Center = center;
581  }
582 
583  itkGetConstMacro(Singular, bool);
584 
585 private:
586  MatrixType m_Matrix; // Matrix of the transformation
587  OutputVectorType m_Offset; // Offset of the transformation
588  mutable InverseMatrixType m_InverseMatrix; // Inverse of the matrix
589  mutable bool m_Singular; // Is m_Inverse singular?
590 
593 
597 }; // class MatrixOffsetTransformBase
598 } // namespace itk
599 
600 #ifndef ITK_MANUAL_INSTANTIATION
601 # include "itkMatrixOffsetTransformBase.hxx"
602 #endif
603 
604 #endif /* itkMatrixOffsetTransformBase_h */
itk::MatrixOffsetTransformBase::GetMatrix
virtual const MatrixType & GetMatrix() const
Definition: itkMatrixOffsetTransformBase.h:246
itk::MatrixOffsetTransformBase::SetVarOffset
void SetVarOffset(const OutputVectorType &offset)
Definition: itkMatrixOffsetTransformBase.h:572
itk::MatrixOffsetTransformBase::m_InverseMatrixMTime
TimeStamp m_InverseMatrixMTime
Definition: itkMatrixOffsetTransformBase.h:596
itk::MatrixOffsetTransformBase::GetVarInverseMatrix
const InverseMatrixType & GetVarInverseMatrix() const
Definition: itkMatrixOffsetTransformBase.h:523
itk::MatrixOffsetTransformBase::InverseMatrixIsOld
bool InverseMatrixIsOld() const
Definition: itkMatrixOffsetTransformBase.h:534
itk::MatrixOffsetTransformBase::SetVarInverseMatrix
void SetVarInverseMatrix(const InverseMatrixType &matrix) const
Definition: itkMatrixOffsetTransformBase.h:528
itk::Transform::InverseJacobianPositionType
vnl_matrix_fixed< ParametersValueType, NInputDimensions, NOutputDimensions > InverseJacobianPositionType
Definition: itkTransform.h:131
itk::Transform::JacobianPositionType
vnl_matrix_fixed< ParametersValueType, NOutputDimensions, NInputDimensions > JacobianPositionType
Definition: itkTransform.h:130
itkMatrix.h
itk::MatrixOffsetTransformBase::m_Singular
bool m_Singular
Definition: itkMatrixOffsetTransformBase.h:589
itk::MatrixOffsetTransformBase::SetVarMatrix
void SetVarMatrix(const MatrixType &matrix)
Definition: itkMatrixOffsetTransformBase.h:553
itk::MatrixOffsetTransformBase< TParametersValueType, NDimensions, NDimensions >::OutputVectorValueType
typename OutputVectorType::ValueType OutputVectorValueType
Definition: itkMatrixOffsetTransformBase.h:150
itk::MatrixOffsetTransformBase::m_InverseMatrix
InverseMatrixType m_InverseMatrix
Definition: itkMatrixOffsetTransformBase.h:588
itk::MatrixOffsetTransformBase::GetTranslation
const OutputVectorType & GetTranslation() const
Definition: itkMatrixOffsetTransformBase.h:347
itk::MatrixOffsetTransformBase::m_Matrix
MatrixType m_Matrix
Definition: itkMatrixOffsetTransformBase.h:583
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itk::MatrixOffsetTransformBase
Matrix and Offset transformation of a vector space (e.g. space coordinates)
Definition: itkMatrixOffsetTransformBase.h:106
itk::DiffusionTensor3D
Represent a diffusion tensor as used in DTI images.
Definition: itkDiffusionTensor3D.h:79
itk::MatrixOffsetTransformBase::SetOffset
void SetOffset(const OutputVectorType &offset)
Definition: itkMatrixOffsetTransformBase.h:260
itk::SmartPointer< Self >
itk::MatrixOffsetTransformBase::SetVarTranslation
void SetVarTranslation(const OutputVectorType &translation)
Definition: itkMatrixOffsetTransformBase.h:563
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::Transform::FixedParametersValueType
typename Superclass::FixedParametersValueType FixedParametersValueType
Definition: itkTransform.h:120
itk::MatrixOffsetTransformBase::SetTranslation
void SetTranslation(const OutputVectorType &translation)
Definition: itkMatrixOffsetTransformBase.h:331
itk::MatrixOffsetTransformBase::m_Center
InputPointType m_Center
Definition: itkMatrixOffsetTransformBase.h:591
itk::MatrixOffsetTransformBase< TParametersValueType, NDimensions, NDimensions >::InputPointValueType
typename InputPointType::ValueType InputPointValueType
Definition: itkMatrixOffsetTransformBase.h:175
itk::MatrixOffsetTransformBase< TParametersValueType, NDimensions, NDimensions >::TranslationValueType
typename TranslationType::ValueType TranslationValueType
Definition: itkMatrixOffsetTransformBase.h:193
itk::SymmetricSecondRankTensor
Represent a symmetric tensor of second rank.
Definition: itkSymmetricSecondRankTensor.h:75
itk::MatrixOffsetTransformBase::m_Offset
OutputVectorType m_Offset
Definition: itkMatrixOffsetTransformBase.h:587
itk::MatrixOffsetTransformBase::GetOffset
const OutputVectorType & GetOffset() const
Definition: itkMatrixOffsetTransformBase.h:275
itk::MatrixOrthogonalityTolerance< float >::GetTolerance
static float GetTolerance()
Definition: itkMatrixOffsetTransformBase.h:57
itk::LightObject
Light weight base class for most itk classes.
Definition: itkLightObject.h:59
itk::Transform::ParametersValueType
typename Superclass::ParametersValueType ParametersValueType
Definition: itkTransform.h:122
itkMacro.h
itk::MatrixOffsetTransformBase::SetCenter
void SetCenter(const InputPointType &center)
Definition: itkMatrixOffsetTransformBase.h:303
itk::MatrixOffsetTransformBase::SetVarCenter
void SetVarCenter(const InputPointType &center)
Definition: itkMatrixOffsetTransformBase.h:578
itk::MatrixOffsetTransformBase< TParametersValueType, NDimensions, NDimensions >::OutputPointValueType
typename OutputPointType::ValueType OutputPointValueType
Definition: itkMatrixOffsetTransformBase.h:177
itk::MatrixOffsetTransformBase::GetTransformCategory
TransformCategoryEnum GetTransformCategory() const override
Definition: itkMatrixOffsetTransformBase.h:210
itk::MatrixOffsetTransformBase< TParametersValueType, NDimensions, NDimensions >::OffsetValueType
typename OffsetType::ValueType OffsetValueType
Definition: itkMatrixOffsetTransformBase.h:189
itk::TimeStamp
Generate a unique, increasing time value.
Definition: itkTimeStamp.h:60
itk::Transform::ScalarType
ParametersValueType ScalarType
Definition: itkTransform.h:126
itk::MatrixOffsetTransformBase::SetMatrix
virtual void SetMatrix(const MatrixType &matrix)
Definition: itkMatrixOffsetTransformBase.h:227
itk::VariableLengthVector
Represents an array whose length can be defined at run-time.
Definition: itkConstantBoundaryCondition.h:28
itk::Matrix< TParametersValueType, Self::OutputSpaceDimension, Self::InputSpaceDimension >
itk::Vector::ValueType
T ValueType
Definition: itkVector.h:71
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk::Transform::InputVnlVectorType
vnl_vector_fixed< TParametersValueType, NInputDimensions > InputVnlVectorType
Definition: itkTransform.h:155
itk::Transform::TransformCategoryEnum
typename Superclass::TransformCategoryEnum TransformCategoryEnum
Definition: itkTransform.h:455
itk::MatrixOffsetTransformBase::m_MatrixMTime
TimeStamp m_MatrixMTime
Definition: itkMatrixOffsetTransformBase.h:595
itk::Transform::ParametersType
typename Superclass::ParametersType ParametersType
Definition: itkTransform.h:121
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::MatrixOffsetTransformBase::m_Translation
OutputVectorType m_Translation
Definition: itkMatrixOffsetTransformBase.h:592
itk::MatrixOffsetTransformBase::IsLinear
bool IsLinear() const override
Definition: itkMatrixOffsetTransformBase.h:491
itk::MatrixOffsetTransformBase::GetCenter
const InputPointType & GetCenter() const
Definition: itkMatrixOffsetTransformBase.h:319
itk::MatrixOrthogonalityTolerance< double >::GetTolerance
static double GetTolerance()
Definition: itkMatrixOffsetTransformBase.h:46
itk::MatrixOrthogonalityTolerance
Definition: itkMatrixOffsetTransformBase.h:39
itk::Math::e
static constexpr double e
Definition: itkMath.h:54
itk::Point
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:53
itk::Matrix< TParametersValueType, Self::OutputSpaceDimension, Self::InputSpaceDimension >::ValueType
TParametersValueType ValueType
Definition: itkMatrix.h:58
itk::Point::ValueType
TCoordRep ValueType
Definition: itkPoint.h:63
itk::Transform
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
itkTransform.h
itk::Transform::FixedParametersType
typename Superclass::FixedParametersType FixedParametersType
Definition: itkTransform.h:119
itk::Array2D
Array2D class representing a 2D array with size defined at construction time.
Definition: itkArray2D.h:45
itk::Transform::InverseTransformBasePointer
typename InverseTransformBaseType::Pointer InverseTransformBasePointer
Definition: itkTransform.h:166
itk::Transform::OutputVnlVectorType
vnl_vector_fixed< TParametersValueType, NOutputDimensions > OutputVnlVectorType
Definition: itkTransform.h:156
itk::MatrixOffsetTransformBase< TParametersValueType, NDimensions, NDimensions >::MatrixValueType
typename MatrixType::ValueType MatrixValueType
Definition: itkMatrixOffsetTransformBase.h:181