ITK  4.6.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 <class TScalar, unsigned int NDimensions>
88  public Transform<TScalar, NDimensions, NDimensions>
89 {
90 public:
96 
98  itkTypeMacro( DisplacementFieldTransform, Transform );
99 
101  itkNewMacro( Self );
102 
105 
108 
112 
115 
118 
121 
125 
129 
132 
134  typedef typename Superclass::InputCovariantVectorType
138 
142 
148 
154 
157 
159  itkStaticConstMacro( Dimension, unsigned int, NDimensions );
160 
165 
168 
177 
180  ScalarType,
181  OutputVectorType::Dimension,
182  Dimension>
184 
189  virtual void SetDisplacementField( DisplacementFieldType* field );
190  itkGetModifiableObjectMacro(DisplacementField, DisplacementFieldType );
192 
195  virtual void SetInverseDisplacementField( DisplacementFieldType * inverseDisplacementField );
196  itkGetModifiableObjectMacro(InverseDisplacementField, DisplacementFieldType );
198 
201  virtual void SetInterpolator( InterpolatorType* interpolator );
202  itkGetModifiableObjectMacro( Interpolator, InterpolatorType );
204 
207  virtual void SetInverseInterpolator( InterpolatorType* interpolator );
208  itkGetModifiableObjectMacro(InverseInterpolator, InterpolatorType );
210 
212  itkGetConstReferenceMacro( DisplacementFieldSetTime, ModifiedTimeType );
213 
216  virtual OutputPointType TransformPoint( const InputPointType& thisPoint )
217  const;
218 
222  {
223  itkExceptionMacro( "TransformVector(Vector) unimplemented, use "
224  "TransformVector(Vector,Point)" );
225  }
227 
229  const
230  {
231  itkExceptionMacro( "TransformVector(Vector) unimplemented, use "
232  "TransformVector(Vector,Point)" );
233  }
234 
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
262  {
263  itkExceptionMacro( "TransformCovariantVector(CovariantVector) "
264  "unimplemented, use TransformCovariantVector(CovariantVector,Point)" );
265  }
267 
269  const InputVectorPixelType &) const
270  {
271  itkExceptionMacro( "TransformCovariantVector(CovariantVector) "
272  "unimplemented, use TransformCovariantVector(CovariantVector,Point)" );
273  }
274 
277  virtual void SetParameters(const ParametersType & params)
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 ParametersType & );
304 
327  JacobianType & j) const
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;
350 
355  virtual void ComputeInverseJacobianWithRespectToPosition(const InputPointType & x, JacobianType & j ) const;
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 );
395 
398  bool GetInverse( Self *inverse ) const;
399 
403 
406  {
408  }
409 
411  {
412  return Dimension;
413  }
414 
415 protected:
416 
418  virtual ~DisplacementFieldTransform();
419  void PrintSelf( std::ostream& os, Indent indent ) const;
420 
424 
428 
432 
436 
437 private:
438  DisplacementFieldTransform( const Self & ); // purposely not implemented
439  void operator=( const Self & ); // purposely not implemented
440 
449  virtual void ComputeJacobianWithRespectToPositionInternal(const IndexType & index, JacobianType & jacobian,
450  bool doInverseJacobian) const;
451 
456  virtual void VerifyFixedParametersInformation();
457 
462  virtual void SetFixedParametersFromDisplacementField() const;
463 
464 };
465 
466 } // end namespace itk
467 
468 #ifndef ITK_MANUAL_INSTANTIATION
469 #include "itkDisplacementFieldTransform.hxx"
470 #endif
471 
472 #endif // __itkDisplacementFieldTransform_h
virtual OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType &) const
Definition: itkTransform.h:238
Superclass::NumberOfParametersType NumberOfParametersType
Definition: itkTransform.h:183
Superclass::ParametersValueType ParametersValueType
virtual NumberOfParametersType GetNumberOfLocalParameters(void) const
virtual void SetInverseDisplacementField(DisplacementFieldType *inverseDisplacementField)
DiffusionTensor3D< TScalar > OutputDiffusionTensor3DType
Definition: itkTransform.h:142
Superclass::RegionType RegionType
Definition: itkImage.h:140
Vector< TScalar, NInputDimensions > InputVectorType
Definition: itkTransform.h:126
DisplacementFieldType::SizeType SizeType
vnl_vector_fixed< TScalar, NInputDimensions > InputVnlVectorType
Definition: itkTransform.h:151
Light weight base class for most itk classes.
virtual OutputCovariantVectorType TransformCovariantVector(const InputCovariantVectorType &) const
virtual OutputVnlVectorType TransformVector(const InputVnlVectorType &) const
Vector< TScalar, NOutputDimensions > OutputVectorType
Definition: itkTransform.h:127
virtual void ComputeJacobianWithRespectToParameters(const InputPointType &, JacobianType &j) const
Superclass::InverseTransformBasePointer InverseTransformBasePointer
InverseTransformBaseType::Pointer InverseTransformBasePointer
Definition: itkTransform.h:164
Represent the size (bounds) of a n-dimensional image.
Definition: itkSize.h:52
Superclass::OutputVectorPixelType OutputVectorPixelType
Superclass::InputVectorPixelType InputVectorPixelType
virtual void SetFixedParametersFromDisplacementField() const
Provides local/dense/high-dimensionaltiy transformation via a a displacement field.
Image< OutputVectorType, Dimension > DisplacementFieldType
unsigned long ModifiedTimeType
Definition: itkIntTypes.h:164
Superclass::OutputVnlVectorType OutputVnlVectorType
DisplacementFieldType::Pointer DisplacementFieldPointer
DisplacementFieldType::PointType PointType
virtual void UpdateTransformParameters(const DerivativeType &update, ScalarType factor=1.0)
virtual void ComputeJacobianWithRespectToPositionInternal(const IndexType &index, JacobianType &jacobian, bool doInverseJacobian) const
Point< TScalar, NInputDimensions > InputPointType
Definition: itkTransform.h:155
void operator=(const Self &)
virtual void SetDisplacementField(DisplacementFieldType *field)
bool GetInverse(Self *inverse) const
virtual void ComputeJacobianWithRespectToParameters(const IndexType &, JacobianType &j) const
ImageVectorOptimizerParametersHelper< ScalarType, OutputVectorType::Dimension, Dimension > OptimizerParametersHelperType
DisplacementFieldType::PixelType PixelType
OutputDiffusionTensor3DType TransformDiffusionTensor(const InputDiffusionTensor3DType &) const
virtual void VerifyFixedParametersInformation()
TPixel PixelType
Definition: itkImage.h:89
VariableLengthVector< TScalar > InputVectorPixelType
Definition: itkTransform.h:131
DisplacementFieldType::ConstPointer DisplacementFieldConstPointer
VariableLengthVector< TScalar > OutputVectorPixelType
Definition: itkTransform.h:132
DiffusionTensor3D< TScalar > InputDiffusionTensor3DType
Definition: itkTransform.h:141
DisplacementFieldType::Pointer m_InverseDisplacementField
OutputVectorPixelType TransformDiffusionTensor(const InputVectorPixelType &) const
CovariantVector< TScalar, NOutputDimensions > OutputCovariantVectorType
Definition: itkTransform.h:148
CovariantVector< TScalar, NInputDimensions > InputCovariantVectorType
Definition: itkTransform.h:146
Transform< TScalar, NDimensions, NDimensions > Superclass
IdentifierType NumberOfParametersType
Superclass::InputDiffusionTensor3DType InputDiffusionTensor3DType
Transform points and vectors from an input space to an output space.
Definition: itkTransform.h:82
virtual void GetInverseJacobianOfForwardFieldWithRespectToPosition(const InputPointType &point, JacobianType &jacobian, bool useSVD=false) const
DisplacementFieldType::SpacingType SpacingType
Superclass::InputVnlVectorType InputVnlVectorType
virtual void ComputeInverseJacobianWithRespectToPosition(const InputPointType &x, JacobianType &j) const
Superclass::ParametersValueType ParametersValueType
Definition: itkTransform.h:119
Class to hold and manage parameters of type Image&lt;Vector&lt;...&gt;,...&gt;, used in Transforms, etc.
Class to hold and manage different parameter types used during optimization.
CovariantVector< ScalarType, OutputDiffusionTensor3DType::Dimension > OutputTensorEigenVectorType
virtual InverseTransformBasePointer GetInverseTransform() const
Superclass::TransformCategoryType TransformCategoryType
Superclass::IndexType IndexType
Definition: itkImage.h:122
Superclass::InputCovariantVectorType InputCovariantVectorType
DisplacementFieldType::DirectionType DirectionType
Superclass::NumberOfParametersType NumberOfParametersType
virtual OutputVectorPixelType TransformCovariantVector(const InputVectorPixelType &) const
virtual void ComputeJacobianWithRespectToPosition(const InputPointType &x, JacobianType &j) const
virtual void Modified() const
Superclass::OutputVectorType OutputVectorType
virtual OutputVectorType TransformVector(const InputVectorType &) const
Definition: itkTransform.h:192
Point< TScalar, NOutputDimensions > OutputPointType
Definition: itkTransform.h:156
virtual void SetFixedParameters(const ParametersType &)
DisplacementFieldType::RegionType RegionType
Superclass::OutputCovariantVectorType OutputCovariantVectorType
virtual void SetInterpolator(InterpolatorType *interpolator)
virtual OutputPointType TransformPoint(const InputPointType &thisPoint) const
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Base class for all vector image interpolaters.
CovariantVector< ScalarType, InputDiffusionTensor3DType::Dimension > InputTensorEigenVectorType
virtual void SetInverseInterpolator(InterpolatorType *interpolator)
Superclass::TransformCategoryType TransformCategoryType
Definition: itkTransform.h:431
SizeValueType Size(void) const
Definition: itkArray.h:124
VectorInterpolateImageFunction< DisplacementFieldType, ScalarType > InterpolatorType
vnl_vector_fixed< TScalar, NOutputDimensions > OutputVnlVectorType
Definition: itkTransform.h:152
virtual OutputVectorType TransformVector(const InputVectorType &) const
virtual OutputVectorPixelType TransformVector(const InputVectorPixelType &) const
virtual void SetParameters(const ParametersType &params)
TScalar ScalarType
Definition: itkTransform.h:115
Superclass::OutputDiffusionTensor3DType OutputDiffusionTensor3DType
A templated class holding a n-Dimensional covariant vector.
DisplacementFieldType::IndexType IndexType
void PrintSelf(std::ostream &os, Indent indent) const
Templated n-dimensional image class.
Definition: itkImage.h:75
virtual TransformCategoryType GetTransformCategory() const
virtual OutputDiffusionTensor3DType TransformDiffusionTensor3D(const InputDiffusionTensor3DType &) const
Definition: itkTransform.h:272
DisplacementFieldType::Pointer m_DisplacementField