ITK  6.0.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  * https://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 VInputDimension = 3, unsigned int VOutputDimension = 3>
106 class ITK_TEMPLATE_EXPORT MatrixOffsetTransformBase
107  : public Transform<TParametersValueType, VInputDimension, VOutputDimension>
108 {
109 public:
110  ITK_DISALLOW_COPY_AND_MOVE(MatrixOffsetTransformBase);
111 
115 
118 
120  itkOverrideGetNameOfClassMacro(MatrixOffsetTransformBase);
121 
123  itkNewMacro(Self);
124 
126  static constexpr unsigned int InputSpaceDimension = VInputDimension;
127  static constexpr unsigned int OutputSpaceDimension = VOutputDimension;
128  static constexpr unsigned int ParametersDimension = VOutputDimension * (VInputDimension + 1);
129 
131  using typename Superclass::FixedParametersType;
132  using typename Superclass::FixedParametersValueType;
133  using typename Superclass::ParametersType;
134  using typename Superclass::ParametersValueType;
135 
137  using typename Superclass::JacobianType;
138  using typename Superclass::JacobianPositionType;
139  using typename Superclass::InverseJacobianPositionType;
140 
142  using typename Superclass::TransformCategoryEnum;
143 
145  using typename Superclass::ScalarType;
146 
151 
155 
156  using typename Superclass::InputVectorPixelType;
157  using typename Superclass::OutputVectorPixelType;
158 
160  using typename Superclass::InputDiffusionTensor3DType;
161  using typename Superclass::OutputDiffusionTensor3DType;
162 
164  using typename Superclass::InputSymmetricSecondRankTensorType;
165  using 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  friend class MatrixOffsetTransformBase<TParametersValueType, VOutputDimension, VInputDimension>;
204 
208  virtual void
209  SetIdentity();
210 
214  TransformCategoryEnum
215  GetTransformCategory() const override
216  {
217  return Self::TransformCategoryEnum::Linear;
218  }
219 
231  virtual void
232  SetMatrix(const MatrixType & matrix)
233  {
234  m_Matrix = matrix;
235  this->ComputeOffset();
236  this->ComputeMatrixParameters();
237  m_MatrixMTime.Modified();
238  this->Modified();
239  return;
240  }
250  virtual const MatrixType &
251  GetMatrix() const
252  {
253  return m_Matrix;
254  }
255 
264  void
265  SetOffset(const OutputVectorType & offset)
266  {
267  m_Offset = offset;
268  this->ComputeTranslation();
269  this->Modified();
270  return;
271  }
279  const OutputVectorType &
280  GetOffset() const
281  {
282  return m_Offset;
283  }
284 
307  void
308  SetCenter(const InputPointType & center)
309  {
310  m_Center = center;
311  this->ComputeOffset();
312  this->Modified();
313  return;
314  }
323  const InputPointType &
324  GetCenter() const
325  {
326  return m_Center;
327  }
328 
335  void
336  SetTranslation(const OutputVectorType & translation)
337  {
338  m_Translation = translation;
339  this->ComputeOffset();
340  this->Modified();
341  return;
342  }
351  const OutputVectorType &
353  {
354  return m_Translation;
355  }
356 
361  void
362  SetParameters(const ParametersType & parameters) override;
363 
365  const ParametersType &
366  GetParameters() const override;
367 
369  void
370  SetFixedParameters(const FixedParametersType &) override;
371 
373  const FixedParametersType &
374  GetFixedParameters() const override;
375 
387  void
388  Compose(const Self * other, bool pre = false);
389 
398  OutputPointType
399  TransformPoint(const InputPointType & point) const override;
400 
401  using Superclass::TransformVector;
402 
403  OutputVectorType
404  TransformVector(const InputVectorType & vect) const override;
405 
406  OutputVnlVectorType
407  TransformVector(const InputVnlVectorType & vect) const override;
408 
409  OutputVectorPixelType
410  TransformVector(const InputVectorPixelType & vect) const override;
411 
412  using Superclass::TransformCovariantVector;
413 
414  OutputCovariantVectorType
415  TransformCovariantVector(const InputCovariantVectorType & vec) const override;
416 
417  OutputVectorPixelType
418  TransformCovariantVector(const InputVectorPixelType & vect) const override;
419 
420  using Superclass::TransformDiffusionTensor3D;
421 
422  OutputDiffusionTensor3DType
423  TransformDiffusionTensor3D(const InputDiffusionTensor3DType & tensor) const override;
424 
425  OutputVectorPixelType
426  TransformDiffusionTensor3D(const InputVectorPixelType & tensor) const override;
427 
428  using Superclass::TransformSymmetricSecondRankTensor;
429  OutputSymmetricSecondRankTensorType
430  TransformSymmetricSecondRankTensor(const InputSymmetricSecondRankTensorType & inputTensor) const override;
431 
432  OutputVectorPixelType
433  TransformSymmetricSecondRankTensor(const InputVectorPixelType & inputTensor) const override;
434 
435 
445  void
446  ComputeJacobianWithRespectToParameters(const InputPointType & p, JacobianType & jacobian) const override;
447 
448 
452  void
453  ComputeJacobianWithRespectToPosition(const InputPointType & x, JacobianPositionType & jac) const override;
454  using Superclass::ComputeJacobianWithRespectToPosition;
455 
459  void
460  ComputeInverseJacobianWithRespectToPosition(const InputPointType & x,
461  InverseJacobianPositionType & jac) const override;
462  using Superclass::ComputeInverseJacobianWithRespectToPosition;
463 
482  bool
483  GetInverse(InverseTransformType * inverse) const;
487  InverseTransformBasePointer
488  GetInverseTransform() const override;
489 
495  bool
496  IsLinear() const override
497  {
498  return true;
499  }
500 
501 protected:
504  const InverseMatrixType &
505  GetInverseMatrix() const;
506 
507 protected:
515 #if !defined(ITK_LEGACY_REMOVE)
516  [[deprecated("Removed unused constructor")]] MatrixOffsetTransformBase(const MatrixType & matrix,
517  const OutputVectorType & offset);
518 #endif
519  explicit MatrixOffsetTransformBase(unsigned int paramDims = ParametersDimension);
523  ~MatrixOffsetTransformBase() override = default;
524 
526  void
527  PrintSelf(std::ostream & os, Indent indent) const override;
528 
529  const InverseMatrixType &
531  {
532  return m_InverseMatrix;
533  }
534  void
536  {
537  m_InverseMatrix = matrix;
538  m_InverseMatrixMTime.Modified();
539  }
540  bool
542  {
543  if (m_MatrixMTime != m_InverseMatrixMTime)
544  {
545  return true;
546  }
547  else
548  {
549  return false;
550  }
551  }
552 
553  virtual void
554  ComputeMatrixParameters();
555 
556  virtual void
557  ComputeMatrix();
558 
559  void
560  SetVarMatrix(const MatrixType & matrix)
561  {
562  m_Matrix = matrix;
563  m_MatrixMTime.Modified();
564  }
565 
566  virtual void
567  ComputeTranslation();
568 
569  void
570  SetVarTranslation(const OutputVectorType & translation)
571  {
572  m_Translation = translation;
573  }
574 
575  virtual void
576  ComputeOffset();
577 
578  void
580  {
581  m_Offset = offset;
582  }
583 
584  void
585  SetVarCenter(const InputPointType & center)
586  {
587  m_Center = center;
588  }
589 
590  itkGetConstMacro(Singular, bool);
591 
592 private:
593  MatrixType m_Matrix{ MatrixType::GetIdentity() }; // Matrix of the transformation
594  OutputVectorType m_Offset{}; // Offset of the transformation
595  mutable InverseMatrixType m_InverseMatrix{ InverseMatrixType::GetIdentity() }; // Inverse of the matrix
596  mutable bool m_Singular{ false }; // Is m_Inverse singular?
597 
598  InputPointType m_Center{};
599  OutputVectorType m_Translation{};
600 
602  TimeStamp m_MatrixMTime{};
603  mutable TimeStamp m_InverseMatrixMTime{};
604 }; // class MatrixOffsetTransformBase
605 } // namespace itk
606 
607 #ifndef ITK_MANUAL_INSTANTIATION
608 # include "itkMatrixOffsetTransformBase.hxx"
609 #endif
610 
611 #endif /* itkMatrixOffsetTransformBase_h */
TransformBaseTemplate
Definition: itkTransformBase.h:85
Pointer
SmartPointer< Self > Pointer
Definition: itkAddImageFilter.h:93
itk::MatrixOffsetTransformBase::GetMatrix
virtual const MatrixType & GetMatrix() const
Definition: itkMatrixOffsetTransformBase.h:251
itk::MatrixOffsetTransformBase::SetVarTranslation
void SetVarTranslation(const OutputVectorType &translation)
Definition: itkMatrixOffsetTransformBase.h:570
itk::MatrixOffsetTransformBase::GetTransformCategory
TransformCategoryEnum GetTransformCategory() const override
Definition: itkMatrixOffsetTransformBase.h:215
itk::MatrixOffsetTransformBase::SetVarCenter
void SetVarCenter(const InputPointType &center)
Definition: itkMatrixOffsetTransformBase.h:585
itk::MatrixOffsetTransformBase::InverseMatrixIsOld
bool InverseMatrixIsOld() const
Definition: itkMatrixOffsetTransformBase.h:541
itk::MatrixOffsetTransformBase::IsLinear
bool IsLinear() const override
Definition: itkMatrixOffsetTransformBase.h:496
itkMatrix.h
itk::Transform::OutputVnlVectorType
vnl_vector_fixed< TParametersValueType, VOutputDimension > OutputVnlVectorType
Definition: itkTransform.h:157
itk::Transform::InverseTransformBasePointer
typename InverseTransformBaseType::Pointer InverseTransformBasePointer
Definition: itkTransform.h:167
itk::MatrixOffsetTransformBase::SetCenter
void SetCenter(const InputPointType &center)
Definition: itkMatrixOffsetTransformBase.h:308
itk::MatrixOffsetTransformBase< TParametersValueType, VDimension, VDimension >::OutputPointValueType
typename OutputPointType::ValueType OutputPointValueType
Definition: itkMatrixOffsetTransformBase.h:177
itk::Vector
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
itk::MatrixOffsetTransformBase< TParametersValueType, VDimension, VDimension >::MatrixValueType
typename MatrixType::ValueType MatrixValueType
Definition: itkMatrixOffsetTransformBase.h:181
itk::MatrixOffsetTransformBase
Matrix and Offset transformation of a vector space (e.g. space coordinates)
Definition: itkMatrixOffsetTransformBase.h:106
itk::Point::ValueType
TCoordRep ValueType
Definition: itkPoint.h:63
itk::SmartPointer< Self >
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::MatrixOffsetTransformBase::SetVarOffset
void SetVarOffset(const OutputVectorType &offset)
Definition: itkMatrixOffsetTransformBase.h:579
itk::MatrixOffsetTransformBase::GetCenter
const InputPointType & GetCenter() const
Definition: itkMatrixOffsetTransformBase.h:324
itk::MatrixOffsetTransformBase::GetVarInverseMatrix
const InverseMatrixType & GetVarInverseMatrix() const
Definition: itkMatrixOffsetTransformBase.h:530
itk::Matrix< TParametersValueType, Self::OutputSpaceDimension, Self::InputSpaceDimension >::ValueType
TParametersValueType ValueType
Definition: itkMatrix.h:59
itk::Vector::ValueType
T ValueType
Definition: itkVector.h:71
itk::MatrixOrthogonalityTolerance< float >::GetTolerance
static float GetTolerance()
Definition: itkMatrixOffsetTransformBase.h:57
itk::MatrixOffsetTransformBase::SetMatrix
virtual void SetMatrix(const MatrixType &matrix)
Definition: itkMatrixOffsetTransformBase.h:232
itk::MatrixOffsetTransformBase< TParametersValueType, VDimension, VDimension >::OutputVectorValueType
typename OutputVectorType::ValueType OutputVectorValueType
Definition: itkMatrixOffsetTransformBase.h:150
itkMacro.h
itk::point
*par Constraints *The filter requires an image with at least two dimensions and a vector *length of at least The theory supports extension to scalar but *the implementation of the itk vector classes do not **The template parameter TRealType must be floating point(float or double) or *a user-defined "real" numerical type with arithmetic operations defined *sufficient to compute derivatives. **\par Performance *This filter will automatically multithread if run with *SetUsePrincipleComponents
itk::MatrixOffsetTransformBase::SetVarInverseMatrix
void SetVarInverseMatrix(const InverseMatrixType &matrix) const
Definition: itkMatrixOffsetTransformBase.h:535
itk::MatrixOffsetTransformBase< TParametersValueType, VDimension, VDimension >::TranslationValueType
typename TranslationType::ValueType TranslationValueType
Definition: itkMatrixOffsetTransformBase.h:193
itk::TimeStamp
Generate a unique, increasing time value.
Definition: itkTimeStamp.h:60
itk::MatrixOffsetTransformBase::GetOffset
const OutputVectorType & GetOffset() const
Definition: itkMatrixOffsetTransformBase.h:280
itk::MatrixOffsetTransformBase< TParametersValueType, VDimension, VDimension >::OffsetValueType
typename OffsetType::ValueType OffsetValueType
Definition: itkMatrixOffsetTransformBase.h:189
itk::MatrixOffsetTransformBase< TParametersValueType, VDimension, VDimension >::InputPointValueType
typename InputPointType::ValueType InputPointValueType
Definition: itkMatrixOffsetTransformBase.h:175
itk::Matrix< TParametersValueType, Self::OutputSpaceDimension, Self::InputSpaceDimension >
itk::MatrixOffsetTransformBase::SetOffset
void SetOffset(const OutputVectorType &offset)
Definition: itkMatrixOffsetTransformBase.h:265
itk::CovariantVector
A templated class holding a n-Dimensional covariant vector.
Definition: itkCovariantVector.h:70
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::MatrixOffsetTransformBase::SetVarMatrix
void SetVarMatrix(const MatrixType &matrix)
Definition: itkMatrixOffsetTransformBase.h:560
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:56
itk::Point
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:53
itk::Transform
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:83
AddImageFilter
Definition: itkAddImageFilter.h:81
itkTransform.h
itk::MatrixOffsetTransformBase::GetTranslation
const OutputVectorType & GetTranslation() const
Definition: itkMatrixOffsetTransformBase.h:352
itk::MatrixOffsetTransformBase::SetTranslation
void SetTranslation(const OutputVectorType &translation)
Definition: itkMatrixOffsetTransformBase.h:336
itk::Transform::InputVnlVectorType
vnl_vector_fixed< TParametersValueType, VInputDimension > InputVnlVectorType
Definition: itkTransform.h:156