ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkDisplacementFieldTransform.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 itkDisplacementFieldTransform_h
19 #define itkDisplacementFieldTransform_h
20 
21 #include "itkTransform.h"
22 
23 #include "itkImage.h"
27 
28 namespace itk
29 {
30 
85 template
86 <typename TParametersValueType, unsigned int NDimensions>
88  public Transform<TParametersValueType, NDimensions, NDimensions>
89 {
90 public:
96 
98  itkTypeMacro( DisplacementFieldTransform, Transform );
99 
101  itkNewMacro( Self );
102 
105 
108 
114 
117 
120 
123 
127 
131 
134 
136  typedef typename Superclass::InputCovariantVectorType
140 
144 
150 
156 
159 
161  itkStaticConstMacro( Dimension, unsigned int, NDimensions );
162 
167 
170 
179 
182  ScalarType,
183  OutputVectorType::Dimension,
184  Dimension>
186 
191  virtual void SetDisplacementField( DisplacementFieldType* field );
192  itkGetModifiableObjectMacro(DisplacementField, DisplacementFieldType );
194 
197  virtual void SetInverseDisplacementField( DisplacementFieldType * inverseDisplacementField );
198  itkGetModifiableObjectMacro(InverseDisplacementField, DisplacementFieldType );
200 
203  virtual void SetInterpolator( InterpolatorType* interpolator );
204  itkGetModifiableObjectMacro( Interpolator, InterpolatorType );
206 
209  virtual void SetInverseInterpolator( InterpolatorType* interpolator );
210  itkGetModifiableObjectMacro(InverseInterpolator, InterpolatorType );
212 
214  itkGetConstReferenceMacro( DisplacementFieldSetTime, ModifiedTimeType );
215 
218  virtual OutputPointType TransformPoint( const InputPointType& thisPoint ) const ITK_OVERRIDE;
219 
222  virtual OutputVectorType TransformVector(const InputVectorType &) const ITK_OVERRIDE
223  {
224  itkExceptionMacro( "TransformVector(Vector) unimplemented, use "
225  "TransformVector(Vector,Point)" );
226  }
228 
229  virtual OutputVectorPixelType TransformVector(const InputVectorPixelType &) const ITK_OVERRIDE
230  {
231  itkExceptionMacro( "TransformVector(Vector) unimplemented, use "
232  "TransformVector(Vector,Point)" );
233  }
234 
235  virtual OutputVnlVectorType TransformVector(const InputVnlVectorType &) const ITK_OVERRIDE
236  {
237  itkExceptionMacro( "TransformVector(Vector) unimplemented, use "
238  "TransformVector(Vector,Point)" );
239  }
240 
244  const InputDiffusionTensor3DType & ) const
245  {
246  itkExceptionMacro( "TransformDiffusionTensor(Tensor) unimplemented, use "
247  "TransformDiffusionTensor(Tensor,Point)" );
248  }
250 
252  const
253  {
254  itkExceptionMacro( "TransformDiffusionTensor(Tensor) unimplemented, use "
255  "TransformDiffusionTensor(Tensor,Point)" );
256  }
257 
261  const InputCovariantVectorType &) const ITK_OVERRIDE
262  {
263  itkExceptionMacro( "TransformCovariantVector(CovariantVector) "
264  "unimplemented, use TransformCovariantVector(CovariantVector,Point)" );
265  }
267 
269  const InputVectorPixelType &) const ITK_OVERRIDE
270  {
271  itkExceptionMacro( "TransformCovariantVector(CovariantVector) "
272  "unimplemented, use TransformCovariantVector(CovariantVector,Point)" );
273  }
274 
277  virtual void SetParameters(const ParametersType & params) ITK_OVERRIDE
278  {
279  if( &(this->m_Parameters) != &params )
280  {
281  if( params.Size() != this->m_Parameters.Size() )
282  {
283  itkExceptionMacro("Input parameters size (" << params.Size()
284  << ") does not match internal size ("
285  << this->m_Parameters.Size() << ").");
286  }
287  /* copy into existing object */
288  this->m_Parameters = params;
289  this->Modified();
290  }
291  }
293 
303  virtual void SetFixedParameters( const FixedParametersType & ) ITK_OVERRIDE;
304 
327  JacobianType & j) const ITK_OVERRIDE
328  {
329  j = this->m_IdentityJacobian;
330  }
332 
340  JacobianType & j) const
341  {
342  j = this->m_IdentityJacobian;
343  }
344 
349  virtual void ComputeJacobianWithRespectToPosition(const InputPointType & x, JacobianType & j ) const ITK_OVERRIDE;
350 
355  virtual void ComputeInverseJacobianWithRespectToPosition(const InputPointType & x, JacobianType & j ) const ITK_OVERRIDE;
356 
361  virtual void ComputeJacobianWithRespectToPosition(const IndexType & x, JacobianType & j ) const;
362 
375  JacobianType & jacobian,
376  bool useSVD = false )
377  const;
378 
391  bool useSVD = false )
392  const;
393 
394  virtual void UpdateTransformParameters( const DerivativeType & update, ScalarType factor = 1.0 ) ITK_OVERRIDE;
395 
398  bool GetInverse( Self *inverse ) const;
399 
402  virtual InverseTransformBasePointer GetInverseTransform() const ITK_OVERRIDE;
403 
404  virtual void SetIdentity();
405 
407  virtual TransformCategoryType GetTransformCategory() const ITK_OVERRIDE
408  {
410  }
411 
412  virtual NumberOfParametersType GetNumberOfLocalParameters(void) const ITK_OVERRIDE
413  {
414  return Dimension;
415  }
416 
424  itkSetMacro(CoordinateTolerance,double);
425  itkGetConstMacro(CoordinateTolerance,double);
427 
435  itkSetMacro(DirectionTolerance,double);
436  itkGetConstMacro(DirectionTolerance,double);
438 
439 protected:
440 
442  virtual ~DisplacementFieldTransform();
443  void PrintSelf( std::ostream& os, Indent indent ) const ITK_OVERRIDE;
444 
448 
452 
456 
460 
461 private:
462  DisplacementFieldTransform( const Self & ); // purposely not implemented
463  void operator=( const Self & ); // purposely not implemented
464 
473  virtual void ComputeJacobianWithRespectToPositionInternal(const IndexType & index, JacobianType & jacobian,
474  bool doInverseJacobian) const;
475 
480  virtual void VerifyFixedParametersInformation();
481 
486  virtual void SetFixedParametersFromDisplacementField() const;
487 
494 
495 };
496 
497 } // end namespace itk
498 
499 #ifndef ITK_MANUAL_INSTANTIATION
500 #include "itkDisplacementFieldTransform.hxx"
501 #endif
502 
503 #endif // itkDisplacementFieldTransform_h
CovariantVector< ScalarType, OutputDiffusionTensor3DType::Dimension > OutputTensorEigenVectorType
Superclass::RegionType RegionType
Definition: itkImage.h:140
Superclass::OutputVectorPixelType OutputVectorPixelType
DisplacementFieldType::IndexType IndexType
Transform< TParametersValueType, NDimensions, NDimensions > Superclass
Light weight base class for most itk classes.
Image< OutputVectorType, Dimension > DisplacementFieldType
DiffusionTensor3D< TParametersValueType > OutputDiffusionTensor3DType
Definition: itkTransform.h:144
Point< TParametersValueType, NInputDimensions > InputPointType
Definition: itkTransform.h:157
IdentifierType NumberOfParametersType
virtual OutputVectorType TransformVector(const InputVectorType &) const override
VectorInterpolateImageFunction< DisplacementFieldType, ScalarType > InterpolatorType
DisplacementFieldType::SizeType SizeType
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
Superclass::NumberOfParametersType NumberOfParametersType
Provides local/dense/high-dimensionaltiy transformation via a a displacement field.
virtual OutputVectorType TransformVector(const InputVectorType &) const
Definition: itkTransform.h:193
virtual void ComputeJacobianWithRespectToParameters(const IndexType &, JacobianType &j) const
DisplacementFieldType::Pointer m_InverseDisplacementField
DisplacementFieldType::Pointer m_DisplacementField
unsigned long ModifiedTimeType
Definition: itkIntTypes.h:164
CovariantVector< TParametersValueType, NOutputDimensions > OutputCovariantVectorType
Definition: itkTransform.h:150
vnl_vector_fixed< TParametersValueType, NOutputDimensions > OutputVnlVectorType
Definition: itkTransform.h:154
virtual TransformCategoryType GetTransformCategory() const override
VariableLengthVector< TParametersValueType > InputVectorPixelType
Definition: itkTransform.h:133
virtual OutputDiffusionTensor3DType TransformDiffusionTensor3D(const InputDiffusionTensor3DType &) const
Definition: itkTransform.h:273
DisplacementFieldType::PointType PointType
Point< TParametersValueType, NOutputDimensions > OutputPointType
Definition: itkTransform.h:158
Superclass::ParametersValueType ParametersValueType
Superclass::InputVectorPixelType InputVectorPixelType
virtual void SetInverseDisplacementField(DisplacementFieldType *inverseDisplacementField)
VariableLengthVector< TParametersValueType > OutputVectorPixelType
Definition: itkTransform.h:134
OutputDiffusionTensor3DType TransformDiffusionTensor(const InputDiffusionTensor3DType &) const
virtual void SetFixedParameters(const FixedParametersType &) override
TPixel PixelType
Definition: itkImage.h:89
virtual void ComputeJacobianWithRespectToParameters(const InputPointType &, JacobianType &j) const override
void operator=(const Self &)
virtual OutputPointType TransformPoint(const InputPointType &thisPoint) const override
Superclass::OutputCovariantVectorType OutputCovariantVectorType
virtual void SetInverseInterpolator(InterpolatorType *interpolator)
virtual void VerifyFixedParametersInformation()
Superclass::InputCovariantVectorType InputCovariantVectorType
ImageVectorOptimizerParametersHelper< ScalarType, OutputVectorType::Dimension, Dimension > OptimizerParametersHelperType
virtual void SetFixedParametersFromDisplacementField() const
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
virtual OutputVectorPixelType TransformVector(const InputVectorPixelType &) const override
Class to hold and manage parameters of type Image&lt;Vector&lt;...&gt;,...&gt;, used in Transforms, etc.
virtual NumberOfParametersType GetNumberOfLocalParameters(void) const override
Class to hold and manage different parameter types used during optimization.
virtual void ComputeJacobianWithRespectToPositionInternal(const IndexType &index, JacobianType &jacobian, bool doInverseJacobian) const
Superclass::IndexType IndexType
Definition: itkImage.h:122
void PrintSelf(std::ostream &os, Indent indent) const override
virtual OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType &) const override
DisplacementFieldType::RegionType RegionType
Superclass::FixedParametersType FixedParametersType
DiffusionTensor3D< TParametersValueType > InputDiffusionTensor3DType
Definition: itkTransform.h:143
virtual OutputVnlVectorType TransformVector(const InputVnlVectorType &) const override
virtual InverseTransformBasePointer GetInverseTransform() const override
DisplacementFieldType::DirectionType DirectionType
vnl_vector_fixed< TParametersValueType, NInputDimensions > InputVnlVectorType
Definition: itkTransform.h:153
DisplacementFieldType::PixelType PixelType
virtual void Modified() const
DisplacementFieldType::SpacingType SpacingType
virtual void GetInverseJacobianOfForwardFieldWithRespectToPosition(const InputPointType &point, JacobianType &jacobian, bool useSVD=false) const
virtual void SetParameters(const ParametersType &params) override
Superclass::FixedParametersValueType FixedParametersValueType
Superclass::OutputVnlVectorType OutputVnlVectorType
bool GetInverse(Self *inverse) const
Superclass::InverseTransformBasePointer InverseTransformBasePointer
virtual void ComputeJacobianWithRespectToPosition(const InputPointType &x, JacobianType &j) const override
ParametersValueType ScalarType
Definition: itkTransform.h:122
Vector< TParametersValueType, NInputDimensions > InputVectorType
Definition: itkTransform.h:128
Superclass::FixedParametersValueType FixedParametersValueType
Definition: itkTransform.h:116
Superclass::OutputVectorType OutputVectorType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
virtual void SetDisplacementField(DisplacementFieldType *field)
virtual OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType &) const
Definition: itkTransform.h:239
Base class for all vector image interpolaters.
Superclass::InputVnlVectorType InputVnlVectorType
SizeValueType Size(void) const
Definition: itkArray.h:124
Vector< TParametersValueType, NOutputDimensions > OutputVectorType
Definition: itkTransform.h:129
DisplacementFieldType::ConstPointer DisplacementFieldConstPointer
Superclass::OutputDiffusionTensor3DType OutputDiffusionTensor3DType
Superclass::InputDiffusionTensor3DType InputDiffusionTensor3DType
DisplacementFieldType::Pointer DisplacementFieldPointer
OutputVectorPixelType TransformDiffusionTensor(const InputVectorPixelType &) const
CovariantVector< ScalarType, InputDiffusionTensor3DType::Dimension > InputTensorEigenVectorType
CovariantVector< TParametersValueType, NInputDimensions > InputCovariantVectorType
Definition: itkTransform.h:148
virtual void ComputeInverseJacobianWithRespectToPosition(const InputPointType &x, JacobianType &j) const override
virtual void SetInterpolator(InterpolatorType *interpolator)
A templated class holding a n-Dimensional covariant vector.
virtual OutputVectorPixelType TransformCovariantVector(const InputVectorPixelType &) const override
virtual void UpdateTransformParameters(const DerivativeType &update, ScalarType factor=1.0) override
Superclass::TransformCategoryType TransformCategoryType
Templated n-dimensional image class.
Definition: itkImage.h:75