Main Page   Groups   Namespace List   Class Hierarchy   Alphabetical List   Compound List   File List   Namespace Members   Compound Members   File Members   Concepts

itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage > Class Template Reference

FEM Image registration filter. More...

#include <itkFEMRegistrationFilter.h>

Inheritance diagram for itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >:

Inheritance graph
[legend]
Collaboration diagram for itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >:

Collaboration graph
[legend]
List of all members.

[NOHEADER]

typedef FieldType::Pointer FieldPointer
FieldPointer ExpandVectorField (ExpandFactorsType *expandFactors, FieldType *f)

Public Types

typedef FEMRegistrationFilter Self
typedef ImageToImageFilter<
TMovingImage, TFixedImage > 
Superclass
typedef SmartPointer< SelfPointer
typedef SmartPointer< const
Self
ConstPointer
typedef TMovingImage MovingImageType
typedef TFixedImage FixedImageType
typedef FixedImageType::PixelType PixelType
typedef FixedImageType::SizeType ImageSizeType
typedef Image< float, itkGetStaticConstMacro(ImageDimension) FloatImageType )
typedef LinearSystemWrapperItpack LinearSystemSolverType
typedef SolverCrankNicolson SolverType
typedef double Float
typedef Load::ArrayType LoadArray
typedef std::vector< typename
LoadLandmark::Pointer
LandmarkArrayType
typedef itk::Vector< float,
itkGetStaticConstMacro(ImageDimension) 
VectorType )
typedef itk::Image< VectorType,
itkGetStaticConstMacro(ImageDimension) 
FieldType )
typedef itk::WarpImageFilter<
MovingImageType, FixedImageType,
FieldType
WarperType
typedef MaterialLinearElasticity MaterialType
typedef itk::ImageRegionIteratorWithIndex<
FixedImageType
ImageIterator
typedef itk::ImageRegionIteratorWithIndex<
FloatImageType
FloatImageIterator
typedef itk::ImageRegionIteratorWithIndex<
FieldType
FieldIterator
typedef itk::VectorIndexSelectionCastImageFilter<
FieldType, FloatImageType
IndexSelectCasterType
typedef double CoordRepType
typedef VectorInterpolateImageFunction<
FieldType, CoordRepType
InterpolatorType
typedef InterpolatorType::Pointer InterpolatorPointer
typedef VectorLinearInterpolateImageFunction<
FieldType, CoordRepType
DefaultInterpolatorType
typedef itk::VectorExpandImageFilter<
FieldType, FieldType
ExpanderType
typedef ExpanderType::ExpandFactorsType ExpandFactorsType
typedef itk::RecursiveMultiResolutionPyramidImageFilter<
FixedImageType, FixedImageType
FixedPyramidType
typedef FiniteDifferenceFunctionLoad<
MovingImageType, FixedImageType
ImageMetricLoadType
typedef PDEDeformableRegistrationFunction<
FixedImageType, MovingImageType,
FieldType
MetricBaseType
typedef MetricBaseType::Pointer MetricBaseTypePointer
enum  Sign { positive = 1, negative = -1 }

Public Member Functions

virtual const char * GetClassName () const
 itkStaticConstMacro (ImageDimension, unsigned int, FixedImageType::ImageDimension)
virtual void SetInterpolator (InterpolatorType *_arg)
virtual InterpolatorTypeGetInterpolator ()
bool ReadConfigFile (const char *)
void RunRegistration ()
void WriteWarpedImage (const char *fn)
void IterativeSolve (SolverType &S)
void MultiResSolve ()
void WarpImage (const MovingImageType *R)
int WriteDisplacementField (unsigned int index)
int WriteDisplacementFieldMultiComponent ()
void SetMovingFile (const char *r)
std::string GetMovingFile ()
void SetFixedFile (const char *t)
std::string GetFixedFile ()
void SetMovingImage (MovingImageType *R)
void SetFixedImage (FixedImageType *T)
MovingImageTypeGetMovingImage ()
MovingImageTypeGetOriginalMovingImage ()
FixedImageTypeGetFixedImage ()
FixedImageTypeGetWarpedImage ()
void ComputeJacobian (float sign=1.0, FieldType *field=NULL, float smooth=0.0)
FloatImageTypeGetJacobianImage ()
FieldTypeGetDeformationField ()
void SetLandmarkFile (const char *l)
void UseLandmarks (bool b)
void EnforceDiffeomorphism (float thresh, SolverType &S, bool onlywriteimages)
void SetResultsFile (const char *r)
void SetResultsFileName (const char *f)
std::string GetResultsFileName ()
void SetDisplacementsFile (const char *r)
void SetNumberOfIntegrationPoints (unsigned int i, unsigned int which=0)
void SetMaximumIterations (unsigned int i, unsigned int which)
void SetTimeStep (Float i)
void SetAlpha (Float a)
void SetEnergyReductionFactor (Float i)
void SetElasticity (Float i, unsigned int which=0)
Float GetElasticity (unsigned int which=0)
void SetRho (Float r, unsigned int which=0)
void SetGamma (Float r, unsigned int which=0)
void SetDescentDirectionMinimize ()
void SetDescentDirectionMaximize ()
void DoLineSearch (unsigned int b)
void DoMultiRes (bool b)
void EmployRegridding (unsigned int b)
void SetLineSearchMaximumIterations (unsigned int f)
void SetWriteDisplacements (bool b)
bool GetWriteDisplacements ()
void SetConfigFileName (const char *f)
std::string GetConfigFileName ()
ImageSizeType GetImageSize ()
void ChooseMetric (float whichmetric)
void SetElement (Element::Pointer e)
void SetMaterial (MaterialType::Pointer m)
void PrintVectorField (unsigned int modnum=1000)
void SetNumLevels (unsigned int i)
void SetMaxLevel (unsigned int i)
void SetTemp (Float i)
void SetDeformationField (FieldType *F)
void SetMeshPixelsPerElementAtEachResolution (unsigned int i, unsigned int which=0)
void SetWidthOfMetricRegion (unsigned int i, unsigned int which=0)
unsigned int GetWidthOfMetricRegion (unsigned int which=0)
MetricBaseTypePointer GetMetric ()
void SetMetric (MetricBaseTypePointer MP)
 FEMRegistrationFilter ()
 ~FEMRegistrationFilter ()

Static Public Member Functions

Pointer New ()

Protected Member Functions

void CreateMesh (double ElementsPerSide, Solver &S, ImageSizeType sz)
void ApplyLoads (SolverType &S, ImageSizeType Isz, double *spacing=NULL)
void ApplyImageLoads (SolverType &S, MovingImageType *i1, FixedImageType *i2)
void CreateLinearSystemSolver ()
Float EvaluateEnergy ()
void InterpolateVectorField (SolverType &S)
FloatImageTypeGetMetricImage (FieldType *F)
void SampleVectorFieldAtNodes (SolverType &S)
Float EvaluateResidual (SolverType &mySolver, Float t)
void FindBracketingTriplet (SolverType &mySolver, Float *a, Float *b, Float *c)
Float GoldenSection (SolverType &mySolver, Float tol=0.01, unsigned int MaxIters=25)
void PrintSelf (std::ostream &os, Indent indent) const
virtual ImageMetricLoadTypeGetLoad ()

Detailed Description

template<class TMovingImage, class TFixedImage>
class itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >

FEM Image registration filter.

The image registration problem is modeled here with the finite element method. Image registration is, in general, an ill-posed problem. Thus, we use an optimization scheme where the optimization criterion is given by a regularized variational energy. The variational energy arises from modeling the image as a physical body on which external forces act. The body is allowed to deform so as to minimize the applied force. The resistance of the physical body to deformation, determined by the physics associated with the body, serves to regularize the solution. The forces applied to the body are, generally, highly non-linear and so the body is allowed to deform slowly and incrementally. The direction it deforms follows the gradient of the potential energy (the force) we define. The potential energies we may choose from are given by the itk image-to-image metrics. The choices and the associated direction of descent are : Mean Squares (minimize), Normalized Cross-Correlation (maximize) Mutual Information (maximize). Note that we have to set the direction (SetDescentDirection) when we choose a metric. The forces driving the problem may also be given by user-supplied landmarks. The corners of the image, in this example, are always pinned. This example is designed for 2D or 3D images. A rectilinear mesh is generated automatically given the correct element type (Quadrilateral or Hexahedral).

Our specific Solver for this example uses trapezoidal time stepping. This is a method for solving a second-order PDE in time. The solution is penalized by the zeroth (mass matrix) and first derivatives (stiffness matrix) of the shape functions. There is an option to perform a line search on the energy after each iteration. Optimal parameter settings require experimentation. The following approach tends to work well : Choose the relative size of density to elasticity (e.g. Rho / E ~= 1.) such that the image deforms locally and slowly. This also affects the stability of the solution. Choose the time step to control the size of the deformation at each step. Choose enough iterations to allow the solution to converge (this may be automated).

Reading images is up to the user. Either set the images using SetMoving/FixedImage or see the ReadImages function.

Note:
This code works for only 2 or 3 dimensions b/c we do not have > 3D elements.

TODO : Keep the full field around (if using re-gridding). Introduce compensation for kinematic non-linearity in time (if using Eulerian frame).

Definition at line 102 of file itkFEMRegistrationFilter.h.


Member Typedef Documentation

template<class TMovingImage, class TFixedImage>
typedef SmartPointer<const Self> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ConstPointer
 

Reimplemented from itk::ImageToImageFilter< TMovingImage, TFixedImage >.

Definition at line 108 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef double itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::CoordRepType
 

Typedef support for the interpolation function Definition at line 145 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef VectorLinearInterpolateImageFunction<FieldType,CoordRepType> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::DefaultInterpolatorType
 

Definition at line 150 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef itk::VectorExpandImageFilter<FieldType,FieldType> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ExpanderType
 

Definition at line 159 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef ExpanderType::ExpandFactorsType itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ExpandFactorsType
 

Definition at line 160 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef itk::ImageRegionIteratorWithIndex<FieldType> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FieldIterator
 

Definition at line 140 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef FieldType::Pointer itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FieldPointer [protected]
 

Re-size the vector field (smaller to larger). Definition at line 440 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef itk::Image<VectorType,itkGetStaticConstMacro(ImageDimension) itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FieldType)
 

Definition at line 135 of file itkFEMRegistrationFilter.h.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetDeformationField().

template<class TMovingImage, class TFixedImage>
typedef TFixedImage itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FixedImageType
 

Definition at line 117 of file itkFEMRegistrationFilter.h.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetFixedImage(), and itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetWarpedImage().

template<class TMovingImage, class TFixedImage>
typedef itk::RecursiveMultiResolutionPyramidImageFilter<FixedImageType,FixedImageType> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FixedPyramidType
 

Definition at line 163 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef double itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::Float
 

Definition at line 129 of file itkFEMRegistrationFilter.h.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetElasticity(), itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetAlpha(), itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetElasticity(), itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetEnergyReductionFactor(), and itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetMaximumIterations().

template<class TMovingImage, class TFixedImage>
typedef itk::ImageRegionIteratorWithIndex<FloatImageType> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FloatImageIterator
 

Definition at line 139 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef Image< float, itkGetStaticConstMacro(ImageDimension) itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FloatImageType)
 

Definition at line 125 of file itkFEMRegistrationFilter.h.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetJacobianImage().

template<class TMovingImage, class TFixedImage>
typedef itk::ImageRegionIteratorWithIndex<FixedImageType> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ImageIterator
 

Definition at line 138 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef FiniteDifferenceFunctionLoad<MovingImageType,FixedImageType> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ImageMetricLoadType
 

Instantiate the load class with the correct image type. Definition at line 171 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef FixedImageType::SizeType itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ImageSizeType
 

Definition at line 119 of file itkFEMRegistrationFilter.h.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetConfigFileName().

template<class TMovingImage, class TFixedImage>
typedef itk::VectorIndexSelectionCastImageFilter<FieldType,FloatImageType> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::IndexSelectCasterType
 

Definition at line 141 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef InterpolatorType::Pointer itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::InterpolatorPointer
 

Definition at line 148 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef VectorInterpolateImageFunction<FieldType,CoordRepType> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::InterpolatorType
 

Definition at line 147 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef std::vector<typename LoadLandmark::Pointer> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::LandmarkArrayType
 

Definition at line 133 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef LinearSystemWrapperItpack itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::LinearSystemSolverType
 

Definition at line 126 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef Load::ArrayType itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::LoadArray
 

Definition at line 130 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef MaterialLinearElasticity itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::MaterialType
 

Definition at line 137 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef PDEDeformableRegistrationFunction<FixedImageType,MovingImageType,FieldType> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::MetricBaseType
 

Definition at line 172 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef MetricBaseType::Pointer itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::MetricBaseTypePointer
 

Definition at line 174 of file itkFEMRegistrationFilter.h.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetImageSize().

template<class TMovingImage, class TFixedImage>
typedef TMovingImage itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::MovingImageType
 

Definition at line 116 of file itkFEMRegistrationFilter.h.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetMovingImage(), and itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetOriginalMovingImage().

template<class TMovingImage, class TFixedImage>
typedef FixedImageType::PixelType itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::PixelType
 

Definition at line 118 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef SmartPointer<Self> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::Pointer
 

Reimplemented from itk::ImageToImageFilter< TMovingImage, TFixedImage >.

Definition at line 107 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef FEMRegistrationFilter itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::Self
 

Standard class typedefs.

Reimplemented from itk::ImageToImageFilter< TMovingImage, TFixedImage >.

Definition at line 105 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef SolverCrankNicolson itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SolverType
 

Definition at line 127 of file itkFEMRegistrationFilter.h.

Referenced by itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::UseLandmarks().

template<class TMovingImage, class TFixedImage>
typedef ImageToImageFilter<TMovingImage, TFixedImage> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::Superclass
 

Reimplemented from itk::ImageToImageFilter< TMovingImage, TFixedImage >.

Definition at line 106 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef itk::Vector<float,itkGetStaticConstMacro(ImageDimension) itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::VectorType)
 

Definition at line 134 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
typedef itk::WarpImageFilter<MovingImageType,FixedImageType, FieldType> itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::WarperType
 

Definition at line 136 of file itkFEMRegistrationFilter.h.


Member Enumeration Documentation

template<class TMovingImage, class TFixedImage>
enum itk::fem::FEMRegistrationFilter::Sign
 

Enumeration values:
positive 
negative 
Definition at line 128 of file itkFEMRegistrationFilter.h.


Constructor & Destructor Documentation

template<class TMovingImage, class TFixedImage>
itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FEMRegistrationFilter  ) 
 

de/constructor

template<class TMovingImage, class TFixedImage>
itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::~FEMRegistrationFilter  ) 
 

de/constructor


Member Function Documentation

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ApplyImageLoads SolverType S,
MovingImageType i1,
FixedImageType i2
[protected]
 

The image loads are entered into the solver.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ApplyLoads SolverType S,
ImageSizeType  Isz,
double *  spacing = NULL
[protected]
 

The non-image loads are entered into the solver.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ChooseMetric float  whichmetric  ) 
 

Choose the metric by parameter : 0= mean squares, 1=cross correlation, 2=pattern intensity, 3 = mutual information.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ComputeJacobian float  sign = 1.0,
FieldType field = NULL,
float  smooth = 0.0
 

Compute the jacobian of the current deformation field.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::CreateLinearSystemSolver  )  [protected]
 

Builds the itpack linear system wrapper with appropriate parameters. Currently undefined

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::CreateMesh double  ElementsPerSide,
Solver S,
ImageSizeType  sz
[protected]
 

This function generates a regular mesh of ElementsPerSide^D size

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::DoLineSearch unsigned int  b  )  [inline]
 

Finds the minimum energy between the current and next solution by linear search. Definition at line 341 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::DoMultiRes bool  b  )  [inline]
 

Sets the use of multi-resolution strategy. The control file always uses multi-res. Definition at line 345 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::EmployRegridding unsigned int  b  )  [inline]
 

Sets the use of multi-resolution strategy. The control file always uses multi-res. Definition at line 348 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::EnforceDiffeomorphism float  thresh,
SolverType S,
bool  onlywriteimages
 

We check the jacobian of the current deformation field. If it is < threshold, we begin diffeomorphism enforcement: 1) Warp the moving image. 2) Set the vector field to zero. 3) Set the warped moving image as the new moving image, resizing if necessary.

template<class TMovingImage, class TFixedImage>
Float itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::EvaluateEnergy  )  [protected]
 

Evaluates the image similarity energy by calling the image metric

template<class TMovingImage, class TFixedImage>
Float itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::EvaluateResidual SolverType mySolver,
Float  t
[protected]
 

template<class TMovingImage, class TFixedImage>
FieldPointer itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ExpandVectorField ExpandFactorsType expandFactors,
FieldType f
[protected]
 

Re-size the vector field (smaller to larger).

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FindBracketingTriplet SolverType mySolver,
Float a,
Float b,
Float c
[protected]
 

template<class TMovingImage, class TFixedImage>
virtual const char* itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetClassName  )  const [virtual]
 

Run-time type information (and related methods)

Reimplemented from itk::ImageToImageFilter< TMovingImage, TFixedImage >.

template<class TMovingImage, class TFixedImage>
std::string itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetConfigFileName  )  [inline]
 

Definition at line 363 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
FieldType* itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetDeformationField  )  [inline]
 

Outputs the FE deformation field interpolated over the entire image domain. Definition at line 241 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FieldType.

template<class TMovingImage, class TFixedImage>
Float itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetElasticity unsigned int  which = 0  )  [inline]
 

Gets the stiffness Matrix weight. Definition at line 326 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::Float.

template<class TMovingImage, class TFixedImage>
std::string itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetFixedFile  )  [inline]
 

Definition at line 212 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
FixedImageType* itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetFixedImage  )  [inline]
 

Definition at line 226 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FixedImageType.

template<class TMovingImage, class TFixedImage>
ImageSizeType itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetImageSize  )  [inline]
 

Definition at line 365 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::MetricBaseTypePointer.

template<class TMovingImage, class TFixedImage>
virtual InterpolatorType* itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetInterpolator  )  [virtual]
 

Get a pointer to the interpolator function.

template<class TMovingImage, class TFixedImage>
FloatImageType* itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetJacobianImage  )  [inline]
 

Get the image that gives the jacobian of the deformation field. Definition at line 237 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FloatImageType.

template<class TMovingImage, class TFixedImage>
virtual ImageMetricLoadType* itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetLoad  )  [protected, virtual]
 

Set the solver's current load.

template<class TMovingImage, class TFixedImage>
MetricBaseTypePointer itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetMetric  )  [inline]
 

Set/Get the Metric. Definition at line 368 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
FloatImageType* itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetMetricImage FieldType F  )  [protected]
 

Calculates the metric over the domain given the vector field.

template<class TMovingImage, class TFixedImage>
std::string itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetMovingFile  )  [inline]
 

Definition at line 208 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
MovingImageType* itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetMovingImage  )  [inline]
 

Definition at line 223 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::MovingImageType.

template<class TMovingImage, class TFixedImage>
MovingImageType* itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetOriginalMovingImage  )  [inline]
 

Definition at line 224 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::MovingImageType.

template<class TMovingImage, class TFixedImage>
std::string itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetResultsFileName  )  [inline]
 

Definition at line 277 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
FixedImageType* itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetWarpedImage  )  [inline]
 

Get the reference image warped to the target image. Must first apply the warp using WarpImage() Definition at line 231 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::FixedImageType.

template<class TMovingImage, class TFixedImage>
unsigned int itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetWidthOfMetricRegion unsigned int  which = 0  )  [inline]
 

The metric region allows one to compute the derivative (force) of the similarity metric using a region of size [i,i] in 2D [i,i,i] in 3D.

Parameters:
i number of elements
which determines the region at a given resolution of the solution process.
Definition at line 302 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
bool itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GetWriteDisplacements  )  [inline]
 

Sets the boolean for writing the displacement field to a file. Definition at line 357 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
Float itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::GoldenSection SolverType mySolver,
Float  tol = 0.01,
unsigned int  MaxIters = 25
[protected]
 

Finds the optimum value between the last two solutions and sets the current solution to that value. Uses Evaluate Residual;

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::InterpolateVectorField SolverType S  )  [protected]
 

Interpolates the vector field over the domain. Our convention is to always keep the vector field at the scale of the original images.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::IterativeSolve SolverType S  ) 
 

The solution loop

template<class TMovingImage, class TFixedImage>
itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::itkStaticConstMacro ImageDimension  ,
unsigned  int,
FixedImageType::ImageDimension 
 

Dimensionality of input and output data is assumed to be the same.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::MultiResSolve  ) 
 

The solution loop for a simple multi-resolution strategy.

template<class TMovingImage, class TFixedImage>
Pointer itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::New  )  [static]
 

Method for creation through the object factory.

Reimplemented from itk::Object.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::PrintSelf std::ostream &  os,
Indent  indent
const [protected, virtual]
 

Methods invoked by Print() to print information about the object including superclasses. Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.

Reimplemented from itk::ImageToImageFilter< TMovingImage, TFixedImage >.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::PrintVectorField unsigned int  modnum = 1000  ) 
 

template<class TMovingImage, class TFixedImage>
bool itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ReadConfigFile const char *   ) 
 

Read the configuration file to set up the example parameters

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::RunRegistration  ) 
 

Call this to register two images.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SampleVectorFieldAtNodes SolverType S  )  [protected]
 

This is used for changing between mesh resolutions.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetAlpha Float  a  )  [inline]
 

Set alpha for the trapezoidal rule (usually 1.0 in our experiments). Definition at line 316 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::Float.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetConfigFileName const char *  f  )  [inline]
 

Sets the file name for the FEM multi-resolution registration. One can also set the parameters in code. Definition at line 361 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::ImageSizeType.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetDeformationField FieldType F  )  [inline]
 

Sets the FE deformation field. Definition at line 244 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetDescentDirectionMaximize  )  [inline]
 

Tries to maximize energy Definition at line 338 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetDescentDirectionMinimize  )  [inline]
 

Tries to minimize energy Definition at line 335 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetDisplacementsFile const char *  r  )  [inline]
 

Sets the filename for the vector field component images. Definition at line 280 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetElasticity Float  i,
unsigned int  which = 0
[inline]
 

Sets the stiffness Matrix weight. Definition at line 323 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::Float.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetElement Element::Pointer  e  )  [inline]
 

This function allows one to set the element and its material externally. Definition at line 377 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetEnergyReductionFactor Float  i  )  [inline]
 

Sets the energy below which we decide the solution has converged. Definition at line 320 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::Float.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetFixedFile const char *  t  )  [inline]
 

Definition at line 210 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetFixedImage FixedImageType T  ) 
 

Define the target (fixed) image.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetGamma Float  r,
unsigned int  which = 0
[inline]
 

Image similarity energy weight Definition at line 332 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
virtual void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetInterpolator InterpolatorType _arg  )  [virtual]
 

Set the interpolator function.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetLandmarkFile const char *  l  )  [inline]
 

These functions control the use of landmark constraints. Currently, landmarks must be read in from a file. Definition at line 253 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetLineSearchMaximumIterations unsigned int  f  )  [inline]
 

This sets the line search's max iterations. Definition at line 351 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetMaterial MaterialType::Pointer  m  )  [inline]
 

This sets the pointer to the material. Definition at line 380 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetMaximumIterations unsigned int  i,
unsigned int  which
[inline]
 

Setting the maximum iterations stops the solution after i iterations regardless of energy.

Parameters:
i number of elements
which determines the resolution of the solution process the call is applied to.
Definition at line 309 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::Float.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetMaxLevel unsigned int  i  )  [inline]
 

Definition at line 385 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetMeshPixelsPerElementAtEachResolution unsigned int  i,
unsigned int  which = 0
[inline]
 

The FEM filter can generate its own mesh for 2 or 3 dimensions, if none is provided. The mesh is generated for quadrilaterals in 2D and hexahedra in 3D. This function sets the number of elements generated along each dimension at the resolution designated by "which". E.g. to generate 10 pixels per element in each dimension in the 1st resolution, use SetMeshResolution(10,0);. Definition at line 288 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetMetric MetricBaseTypePointer  MP  )  [inline]
 

Set/Get the Metric. Definition at line 369 of file itkFEMRegistrationFilter.h.

References itk::fem::Element::Pointer.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetMovingFile const char *  r  )  [inline]
 

One can set the reference file names to read images from files Definition at line 206 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetMovingImage MovingImageType R  ) 
 

Define the reference (moving) image.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetNumberOfIntegrationPoints unsigned int  i,
unsigned int  which = 0
[inline]
 

This determines the number of integration points to use at each resolution. These integration points are used to generate the force. The actual number used will be i^d, where d is the number of parameters in the elements local domain. Definition at line 294 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetNumLevels unsigned int  i  )  [inline]
 

Definition at line 384 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetResultsFile const char *  r  )  [inline]
 

The warped reference image will be written to this file name with the extension "11.img" appended to it. One can also output the image after every iteration, yielding result11.img, result12.img, etc. by uncommenting the code at the end of IterativeSolve. Definition at line 273 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetResultsFileName const char *  f  )  [inline]
 

Definition at line 275 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetRho Float  r,
unsigned int  which = 0
[inline]
 

Mass matrix weight Definition at line 329 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetTemp Float  i  )  [inline]
 

Definition at line 387 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMOF.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetTimeStep Float  i  )  [inline]
 

Setting the time step - usually 1.0. We prefer to use rho to control step sizes. Definition at line 313 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetWidthOfMetricRegion unsigned int  i,
unsigned int  which = 0
[inline]
 

The metric region allows one to compute the derivative (force) of the similarity metric using a region of size [i,i] in 2D [i,i,i] in 3D.

Parameters:
i number of elements
which determines the region at a given resolution of the solution process.
Definition at line 301 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SetWriteDisplacements bool  b  )  [inline]
 

Sets the boolean for writing the displacement field to a file. Definition at line 354 of file itkFEMRegistrationFilter.h.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::UseLandmarks bool  b  )  [inline]
 

This determines if the landmark file will be read Definition at line 256 of file itkFEMRegistrationFilter.h.

References itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::SolverType.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::WarpImage const MovingImageType R  ) 
 

Applies the warp to the input image.

template<class TMovingImage, class TFixedImage>
int itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::WriteDisplacementField unsigned int  index  ) 
 

Writes the displacement field to a file.

template<class TMovingImage, class TFixedImage>
int itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::WriteDisplacementFieldMultiComponent  ) 
 

Writes the displacement field to a file as a single volume with multiple components.

template<class TMovingImage, class TFixedImage>
void itk::fem::FEMRegistrationFilter< TMovingImage, TFixedImage >::WriteWarpedImage const char *  fn  ) 
 

Call this to write out images - a counter is attached to the file name so we can output a numbered sequence tracking the deformation.


The documentation for this class was generated from the following file:
Generated at Sun Apr 1 03:17:17 2007 for ITK by doxygen 1.3.8 written by Dimitri van Heesch, © 1997-2000