ITK  4.13.0
Insight Segmentation and Registration Toolkit
itkMatrixOffsetTransformBase.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 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  * matricies 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 GetTolerance() { return 1e-10; }
46 };
47 
48 template <>
49 class ITK_TEMPLATE_EXPORT MatrixOrthogonalityTolerance<float>
50 {
51 public:
52  static float GetTolerance() { return 1e-5f; }
53 };
54 
96 template<typename TParametersValueType=double,
97  unsigned int NInputDimensions = 3,
98  unsigned int NOutputDimensions = 3>
99 class ITK_TEMPLATE_EXPORT MatrixOffsetTransformBase :
100  public Transform<TParametersValueType, NInputDimensions, NOutputDimensions>
101 {
102 public:
105  typedef Transform<TParametersValueType,
106  NInputDimensions,
107  NOutputDimensions> Superclass;
108 
111 
113  itkTypeMacro(MatrixOffsetTransformBase, Transform);
114 
116  itkNewMacro(Self);
117 
119  itkStaticConstMacro(InputSpaceDimension, unsigned int, NInputDimensions);
120  itkStaticConstMacro(OutputSpaceDimension, unsigned int, NOutputDimensions);
121  itkStaticConstMacro( ParametersDimension, unsigned int,
122  NOutputDimensions * ( NInputDimensions + 1 ) );
124 
126  typedef typename Superclass::FixedParametersType FixedParametersType;
127  typedef typename Superclass::FixedParametersValueType FixedParametersValueType;
128  typedef typename Superclass::ParametersType ParametersType;
129  typedef typename Superclass::ParametersValueType ParametersValueType;
130 
132  typedef typename Superclass::JacobianType JacobianType;
133 
135  typedef typename Superclass::TransformCategoryType TransformCategoryType;
136 
138  typedef typename Superclass::ScalarType ScalarType;
139 
141  typedef Vector<TParametersValueType,
142  itkGetStaticConstMacro(InputSpaceDimension)> InputVectorType;
143  typedef Vector<TParametersValueType,
144  itkGetStaticConstMacro(OutputSpaceDimension)> OutputVectorType;
147 
149  typedef CovariantVector<TParametersValueType,
150  itkGetStaticConstMacro(InputSpaceDimension)>
152  typedef CovariantVector<TParametersValueType,
153  itkGetStaticConstMacro(OutputSpaceDimension)>
156 
157  typedef typename Superclass::InputVectorPixelType InputVectorPixelType;
158  typedef typename Superclass::OutputVectorPixelType OutputVectorPixelType;
159 
161  typedef typename Superclass::InputDiffusionTensor3DType
163  typedef typename Superclass::OutputDiffusionTensor3DType
165 
167  typedef typename Superclass::InputSymmetricSecondRankTensorType
169  typedef typename Superclass::OutputSymmetricSecondRankTensorType
171 
174 
176  typedef vnl_vector_fixed<TParametersValueType,
177  itkGetStaticConstMacro(InputSpaceDimension)>
179  typedef vnl_vector_fixed<TParametersValueType,
180  itkGetStaticConstMacro(OutputSpaceDimension)>
183 
185  typedef Point<TParametersValueType,
186  itkGetStaticConstMacro(InputSpaceDimension)>
189  typedef Point<TParametersValueType,
190  itkGetStaticConstMacro(OutputSpaceDimension)>
194 
196  typedef Matrix<TParametersValueType, itkGetStaticConstMacro(OutputSpaceDimension),
197  itkGetStaticConstMacro(InputSpaceDimension)>
201 
203  typedef Matrix<TParametersValueType, itkGetStaticConstMacro(InputSpaceDimension),
204  itkGetStaticConstMacro(OutputSpaceDimension)>
206 
208 
211 
213 
215 
218  typedef typename Superclass::InverseTransformBaseType InverseTransformBaseType;
220 
224  virtual void SetIdentity();
225 
229  virtual TransformCategoryType GetTransformCategory() const ITK_OVERRIDE
230  {
231  return Self::Linear;
232  }
233 
245  virtual void SetMatrix(const MatrixType & matrix)
246  {
247  m_Matrix = matrix; this->ComputeOffset();
248  this->ComputeMatrixParameters();
249  m_MatrixMTime.Modified(); this->Modified(); return;
250  }
252 
260  virtual const MatrixType & GetMatrix() const
261  {
262  return m_Matrix;
263  }
264 
273  void SetOffset(const OutputVectorType & offset)
274  {
275  m_Offset = offset; this->ComputeTranslation();
276  this->Modified(); return;
277  }
279 
285  const OutputVectorType & GetOffset() const
286  {
287  return m_Offset;
288  }
289 
312  void SetCenter(const InputPointType & center)
313  {
314  m_Center = center; this->ComputeOffset();
315  this->Modified(); return;
316  }
318 
325  const InputPointType & GetCenter() const
326  {
327  return m_Center;
328  }
329 
336  void SetTranslation(const OutputVectorType & translation)
337  {
338  m_Translation = translation; this->ComputeOffset();
339  this->Modified(); return;
340  }
342 
350  {
351  return m_Translation;
352  }
353 
358  void SetParameters(const ParametersType & parameters) ITK_OVERRIDE;
359 
361  const ParametersType & GetParameters() const ITK_OVERRIDE;
362 
364  virtual void SetFixedParameters(const FixedParametersType &) ITK_OVERRIDE;
365 
367  virtual const FixedParametersType & GetFixedParameters() const ITK_OVERRIDE;
368 
380  void Compose(const Self *other, bool pre = 0);
381 
390  OutputPointType TransformPoint(const InputPointType & point) const ITK_OVERRIDE;
391 
392  using Superclass::TransformVector;
393 
394  OutputVectorType TransformVector(const InputVectorType & vector) const ITK_OVERRIDE;
395 
396  OutputVnlVectorType TransformVector(const InputVnlVectorType & vector) const ITK_OVERRIDE;
397 
398  OutputVectorPixelType TransformVector(const InputVectorPixelType & vector) const ITK_OVERRIDE;
399 
400  using Superclass::TransformCovariantVector;
401 
402  OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType & vector) const ITK_OVERRIDE;
403 
404  OutputVectorPixelType TransformCovariantVector(const InputVectorPixelType & vector) const ITK_OVERRIDE;
405 
406  using Superclass::TransformDiffusionTensor3D;
407 
408  OutputDiffusionTensor3DType TransformDiffusionTensor3D(const InputDiffusionTensor3DType & tensor) const ITK_OVERRIDE;
409 
410  OutputVectorPixelType TransformDiffusionTensor3D(const InputVectorPixelType & tensor ) const ITK_OVERRIDE;
411 
412  using Superclass::TransformSymmetricSecondRankTensor;
413  OutputSymmetricSecondRankTensorType TransformSymmetricSecondRankTensor( const InputSymmetricSecondRankTensorType & tensor ) const ITK_OVERRIDE;
414 
415  OutputVectorPixelType TransformSymmetricSecondRankTensor( const InputVectorPixelType & tensor ) const ITK_OVERRIDE;
416 
426  virtual void ComputeJacobianWithRespectToParameters(const InputPointType & x, JacobianType & j) const ITK_OVERRIDE;
427 
431  virtual void ComputeJacobianWithRespectToPosition(const InputPointType & x, JacobianType & jac) const ITK_OVERRIDE;
432 
436  virtual void ComputeInverseJacobianWithRespectToPosition(const InputPointType & x, JacobianType & jac) const ITK_OVERRIDE;
437 
456  bool GetInverse(Self *inverse) const;
458 
460  virtual InverseTransformBasePointer GetInverseTransform() const ITK_OVERRIDE;
461 
467  virtual bool IsLinear() const ITK_OVERRIDE
468  {
469  return true;
470  }
471 
472 #if !defined(ITK_LEGACY_REMOVE)
473 
474 public:
475 #else
476 
477 protected:
478 #endif
479 
481  const InverseMatrixType & GetInverseMatrix() const;
482 
483 protected:
491  MatrixOffsetTransformBase(const MatrixType & matrix, const OutputVectorType & offset);
492  MatrixOffsetTransformBase(unsigned int paramDims);
495 
497  virtual ~MatrixOffsetTransformBase() ITK_OVERRIDE;
498 
500  virtual void PrintSelf(std::ostream & s, Indent indent) const ITK_OVERRIDE;
501 
502  const InverseMatrixType & GetVarInverseMatrix() const
503  {
504  return m_InverseMatrix;
505  }
506  void SetVarInverseMatrix(const InverseMatrixType & matrix) const
507  {
508  m_InverseMatrix = matrix; m_InverseMatrixMTime.Modified();
509  }
510  bool InverseMatrixIsOld() const
511  {
512  if( m_MatrixMTime != m_InverseMatrixMTime )
513  {
514  return true;
515  }
516  else
517  {
518  return false;
519  }
520  }
521 
522  virtual void ComputeMatrixParameters();
523 
524  virtual void ComputeMatrix();
525 
526  void SetVarMatrix(const MatrixType & matrix)
527  {
528  m_Matrix = matrix; m_MatrixMTime.Modified();
529  }
530 
531  virtual void ComputeTranslation();
532 
533  void SetVarTranslation(const OutputVectorType & translation)
534  {
535  m_Translation = translation;
536  }
537 
538  virtual void ComputeOffset();
539 
540  void SetVarOffset(const OutputVectorType & offset)
541  {
542  m_Offset = offset;
543  }
544 
545  void SetVarCenter(const InputPointType & center)
546  {
547  m_Center = center;
548  }
549 
550  itkGetConstMacro(Singular, bool);
551 private:
552 
553  MatrixOffsetTransformBase(const Self & other);
554  const Self & operator=(const Self &);
555 
556  MatrixType m_Matrix; // Matrix of the transformation
557  OutputVectorType m_Offset; // Offset of the transformation
558  mutable InverseMatrixType m_InverseMatrix; // Inverse of the matrix
559  mutable bool m_Singular; // Is m_Inverse singular?
560 
563 
567 }; // class MatrixOffsetTransformBase
568 } // namespace itk
569 
570 #ifndef ITK_MANUAL_INSTANTIATION
571 #include "itkMatrixOffsetTransformBase.hxx"
572 #endif
573 
574 #endif /* itkMatrixOffsetTransformBase_h */
TCoordRep ValueType
Definition: itkPoint.h:61
A templated class holding a M x N size Matrix.
Definition: itkMatrix.h:53
Superclass::InputDiffusionTensor3DType InputDiffusionTensor3DType
Point< TParametersValueType, itkGetStaticConstMacro(OutputSpaceDimension)> OutputPointType
void SetVarCenter(const InputPointType &center)
CovariantVector< TParametersValueType, itkGetStaticConstMacro(InputSpaceDimension)> InputCovariantVectorType
Matrix< TParametersValueType, itkGetStaticConstMacro(OutputSpaceDimension), itkGetStaticConstMacro(InputSpaceDimension)> MatrixType
Light weight base class for most itk classes.
Superclass::OutputVectorPixelType OutputVectorPixelType
void SetCenter(const InputPointType &center)
CovariantVector< TParametersValueType, itkGetStaticConstMacro(OutputSpaceDimension)> OutputCovariantVectorType
Matrix and Offset transformation of a vector space (e.g. space coordinates)
Superclass::FixedParametersType FixedParametersType
Superclass::ParametersValueType ParametersValueType
Superclass::InputSymmetricSecondRankTensorType InputSymmetricSecondRankTensorType
vnl_vector_fixed< TParametersValueType, itkGetStaticConstMacro(OutputSpaceDimension)> OutputVnlVectorType
OutputPointType::ValueType OutputPointValueType
Superclass::TransformCategoryType TransformCategoryType
void SetVarTranslation(const OutputVectorType &translation)
const OutputVectorType & GetTranslation() const
Transform< TParametersValueType, NInputDimensions, NOutputDimensions > Superclass
TranslationType::ValueType TranslationValueType
Superclass::OutputSymmetricSecondRankTensorType OutputSymmetricSecondRankTensorType
void SetVarMatrix(const MatrixType &matrix)
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
Matrix< TParametersValueType, itkGetStaticConstMacro(InputSpaceDimension), itkGetStaticConstMacro(OutputSpaceDimension)> InverseMatrixType
A templated class holding a n-Dimensional vector.
Definition: itkVector.h:62
Superclass::FixedParametersValueType FixedParametersValueType
void SetTranslation(const OutputVectorType &translation)
virtual TransformCategoryType GetTransformCategory() const override
Generate a unique, increasing time value.
Definition: itkTimeStamp.h:59
Vector< TParametersValueType, itkGetStaticConstMacro(OutputSpaceDimension)> OutputVectorType
OutputVectorType::ValueType OutputVectorValueType
Superclass::OutputDiffusionTensor3DType OutputDiffusionTensor3DType
virtual void SetMatrix(const MatrixType &matrix)
virtual const MatrixType & GetMatrix() const
Control indentation during Print() invocation.
Definition: itkIndent.h:49
vnl_vector_fixed< TParametersValueType, itkGetStaticConstMacro(InputSpaceDimension)> InputVnlVectorType
void SetVarInverseMatrix(const InverseMatrixType &matrix) const
const OutputVectorType & GetOffset() const
InverseTransformBaseType::Pointer InverseTransformBasePointer
const InputPointType & GetCenter() const
Point< TParametersValueType, itkGetStaticConstMacro(InputSpaceDimension)> InputPointType
static ITK_CONSTEXPR_VAR double e
The base of the natural logarithm or Euler&#39;s number
Definition: itkMath.h:56
A templated class holding a geometric point in n-Dimensional space.
Definition: itkPoint.h:52
A templated class holding a n-Dimensional covariant vector.
void SetVarOffset(const OutputVectorType &offset)
InputPointType::ValueType InputPointValueType
Superclass::InputVectorPixelType InputVectorPixelType
void SetOffset(const OutputVectorType &offset)
CovariantVector< TParametersValueType, InputDiffusionTensor3DType::Dimension > InputTensorEigenVectorType
Vector< TParametersValueType, itkGetStaticConstMacro(InputSpaceDimension)> InputVectorType
Superclass::InverseTransformBaseType InverseTransformBaseType