SimpleITK  2.0.0
Classes | Public Types | Public Member Functions | Protected Member Functions | Private Types | Private Attributes | Friends | List of all members
itk::simple::ImageRegistrationMethod Class Reference

An interface method to the modular ITKv4 registration framework. More...

#include <sitkImageRegistrationMethod.h>

+ Inheritance diagram for itk::simple::ImageRegistrationMethod:
+ Collaboration diagram for itk::simple::ImageRegistrationMethod:

Classes

struct  EvaluateMemberFunctionAddressor
 

Public Types

enum  EstimateLearningRateType {
  Never,
  Once,
  EachIteration
}
 
enum  MetricSamplingStrategyType {
  NONE,
  REGULAR,
  RANDOM
}
 
using Self = ImageRegistrationMethod
 
using Superclass = ProcessObject
 
- Public Types inherited from itk::simple::ProcessObject
using Self = ProcessObject
 

Public Member Functions

Transform Execute (const Image &fixed, const Image &moving)
 Optimize the configured registration problem. More...
 
unsigned int GetCurrentLevel () const
 
uint64_t GetMetricNumberOfValidPoints () const
 
const std::vector< double > & GetMetricSamplingPercentagePerLevel () const
 Get the percentage of pixels used for metric evaluation.
More...
 
double GetMetricValue () const
 
std::string GetName () const override
 
double GetOptimizerConvergenceValue () const
 
unsigned int GetOptimizerIteration () const
 
double GetOptimizerLearningRate () const
 
std::vector< double > GetOptimizerPosition () const
 
std::vector< double > GetOptimizerScales () const
 Get the OptimizerScales. More...
 
std::string GetOptimizerStopConditionDescription () const
 
 ImageRegistrationMethod ()
 
double MetricEvaluate (const Image &fixed, const Image &moving)
 Get the value of the metric given the state of the method. More...
 
void SetInitialTransformAsBSpline (BSplineTransform &transform, bool inPlace=true, const std::vector< unsigned int > &scaleFactors=std::vector< unsigned int >())
 Set an initial BSpline transform to optimize. More...
 
SelfSetMetricAsANTSNeighborhoodCorrelation (unsigned int radius)
 Use normalized cross correlation using a small neighborhood for each voxel between two images, with speed optimizations for dense registration. More...
 
SelfSetMetricAsCorrelation ()
 Use negative normalized cross correlation image metric. More...
 
SelfSetMetricAsDemons (double intensityDifferenceThreshold=0.001)
 Use demons image metric. More...
 
SelfSetMetricAsJointHistogramMutualInformation (unsigned int numberOfHistogramBins=20, double varianceForJointPDFSmoothing=1.5)
 Use mutual information between two images. More...
 
SelfSetMetricAsMattesMutualInformation (unsigned int numberOfHistogramBins=50)
 Use the mutual information between two images to be registered using the method of Mattes et al. More...
 
SelfSetMetricAsMeanSquares ()
 Use negative means squares image metric. More...
 
SelfSetMetricFixedMask (const Image &binaryMask)
 Set an image mask in order to restrict the sampled points for the metric. More...
 
SelfSetMetricMovingMask (const Image &binaryMask)
 Set an image mask in order to restrict the sampled points for the metric in the moving image space. More...
 
SelfSetMetricSamplingStrategy (MetricSamplingStrategyType strategy)
 Set sampling strategy for sample generation. More...
 
SelfSetOptimizerAsAmoeba (double simplexDelta, unsigned int numberOfIterations, double parametersConvergenceTolerance=1e-8, double functionConvergenceTolerance=1e-4, bool withRestarts=false)
 Set optimizer to Nelder-Mead downhill simplex algorithm. More...
 
SelfSetOptimizerAsConjugateGradientLineSearch (double learningRate, unsigned int numberOfIterations, double convergenceMinimumValue=1e-6, unsigned int convergenceWindowSize=10, double lineSearchLowerLimit=0, double lineSearchUpperLimit=5.0, double lineSearchEpsilon=0.01, unsigned int lineSearchMaximumIterations=20, EstimateLearningRateType estimateLearningRate=Once, double maximumStepSizeInPhysicalUnits=0.0)
 Conjugate gradient descent optimizer with a golden section line search for nonlinear optimization. More...
 
SelfSetOptimizerAsExhaustive (const std::vector< unsigned int > &numberOfSteps, double stepLength=1.0)
 Set the optimizer to sample the metric at regular steps. More...
 
SelfSetOptimizerAsGradientDescent (double learningRate, unsigned int numberOfIterations, double convergenceMinimumValue=1e-6, unsigned int convergenceWindowSize=10, EstimateLearningRateType estimateLearningRate=Once, double maximumStepSizeInPhysicalUnits=0.0)
 Gradient descent optimizer. More...
 
SelfSetOptimizerAsGradientDescentLineSearch (double learningRate, unsigned int numberOfIterations, double convergenceMinimumValue=1e-6, unsigned int convergenceWindowSize=10, double lineSearchLowerLimit=0, double lineSearchUpperLimit=5.0, double lineSearchEpsilon=0.01, unsigned int lineSearchMaximumIterations=20, EstimateLearningRateType estimateLearningRate=Once, double maximumStepSizeInPhysicalUnits=0.0)
 Gradient descent optimizer with a golden section line search. More...
 
SelfSetOptimizerAsLBFGS2 (double solutionAccuracy=1e-5, unsigned int numberOfIterations=0, unsigned int hessianApproximateAccuracy=6, unsigned int deltaConvergenceDistance=0, double deltaConvergenceTolerance=1e-5, unsigned int lineSearchMaximumEvaluations=40, double lineSearchMinimumStep=1e-20, double lineSearchMaximumStep=1e20, double lineSearchAccuracy=1e-4)
 Limited memory Broyden Fletcher Goldfarb Shannon minimization without bounds. More...
 
SelfSetOptimizerAsLBFGSB (double gradientConvergenceTolerance=1e-5, unsigned int numberOfIterations=500, unsigned int maximumNumberOfCorrections=5, unsigned int maximumNumberOfFunctionEvaluations=2000, double costFunctionConvergenceFactor=1e+7, double lowerBound=std::numeric_limits< double >::min(), double upperBound=std::numeric_limits< double >::max(), bool trace=false)
 Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds. More...
 
SelfSetOptimizerAsOnePlusOneEvolutionary (unsigned int numberOfIterations=100, double epsilon=1.5e-4, double initialRadius=1.01, double growthFactor=-1.0, double shrinkFactor=-1.0, unsigned int seed=sitkWallClock)
 1+1 evolutionary optimizer strategy. More...
 
SelfSetOptimizerAsPowell (unsigned int numberOfIterations=100, unsigned int maximumLineIterations=100, double stepLength=1, double stepTolerance=1e-6, double valueTolerance=1e-6)
 Powell optimization using Brent line search. More...
 
SelfSetOptimizerAsRegularStepGradientDescent (double learningRate, double minStep, unsigned int numberOfIterations, double relaxationFactor=0.5, double gradientMagnitudeTolerance=1e-4, EstimateLearningRateType estimateLearningRate=Never, double maximumStepSizeInPhysicalUnits=0.0)
 Regular Step Gradient descent optimizer. More...
 
SelfSetOptimizerScales (const std::vector< double > &scales)
 Manually set per parameter weighting for the transform parameters. More...
 
SelfSetOptimizerScalesFromIndexShift (unsigned int centralRegionRadius=5, double smallParameterVariation=0.01)
 Estimate scales from maximum voxel shift in index space cause by parameter change. More...
 
SelfSetOptimizerScalesFromJacobian (unsigned int centralRegionRadius=5)
 Estimate scales from Jacobian norms. More...
 
SelfSetOptimizerScalesFromPhysicalShift (unsigned int centralRegionRadius=5, double smallParameterVariation=0.01)
 Estimating scales of transform parameters a step sizes, from the maximum voxel shift in physical space caused by a parameter change. More...
 
SelfSetShrinkFactorsPerLevel (const std::vector< unsigned int > &shrinkFactors)
 Set the isotropic shrink factors for each level. More...
 
SelfSetSmoothingSigmasPerLevel (const std::vector< double > &smoothingSigmas)
 Set the sigmas of Gaussian used for smoothing. More...
 
std::string ToString () const override
 Print the information about the object to a string. More...
 
 ~ImageRegistrationMethod () override
 
InterpolatorEnum GetInterpolator ()
 Set and get the interpolator to use. More...
 
SelfSetInterpolator (InterpolatorEnum Interpolator)
 Set and get the interpolator to use. More...
 
SelfSetInitialTransform (const Transform &transform)
 Set the initial transform and parameters to optimize. More...
 
SelfSetInitialTransform (Transform &transform, bool inPlace=true)
 Set the initial transform and parameters to optimize. More...
 
Transform GetInitialTransform ()
 Set the initial transform and parameters to optimize. More...
 
bool GetInitialTransformInPlace () const
 Set the initial transform and parameters to optimize. More...
 
SelfSetMovingInitialTransform (const Transform &transform)
 Set a fixed transform component towards moving domain. More...
 
Transform GetMovingInitialTransform () const
 Set a fixed transform component towards moving domain. More...
 
SelfSetFixedInitialTransform (const Transform &transform)
 Set transform mapping to the fixed domain. More...
 
Transform GetFixedInitialTransform () const
 Set transform mapping to the fixed domain. More...
 
SelfSetVirtualDomain (const std::vector< uint32_t > &virtualSize, const std::vector< double > &virtualOrigin, const std::vector< double > &virtualSpacing, const std::vector< double > &virtualDirection)
 Set the virtual domain used for sampling. More...
 
SelfSetVirtualDomainFromImage (const Image &virtualImage)
 Set the virtual domain used for sampling. More...
 
SelfSetOptimizerWeights (const std::vector< double > &weights)
 A per parameter weighting array for the optimizer. More...
 
std::vector< double > GetOptimizerWeights () const
 A per parameter weighting array for the optimizer. More...
 
SelfSetMetricSamplingPercentage (double percentage, unsigned int seed=sitkWallClock)
 Set percentage of pixels sampled for metric evaluation. More...
 
SelfSetMetricSamplingPercentagePerLevel (const std::vector< double > &percentage, unsigned int seed=sitkWallClock)
 Set percentage of pixels sampled for metric evaluation. More...
 
SelfSetMetricUseFixedImageGradientFilter (bool)
 Enable image gradient computation by a filter. More...
 
SelfMetricUseFixedImageGradientFilterOn ()
 Enable image gradient computation by a filter. More...
 
SelfMetricUseFixedImageGradientFilterOff ()
 Enable image gradient computation by a filter. More...
 
SelfSetMetricUseMovingImageGradientFilter (bool)
 Enable image gradient computation by a filter. More...
 
SelfMetricUseMovingImageGradientFilterOn ()
 Enable image gradient computation by a filter. More...
 
SelfMetricUseMovingImageGradientFilterOff ()
 Enable image gradient computation by a filter. More...
 
SelfSetSmoothingSigmasAreSpecifiedInPhysicalUnits (bool arg)
 Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels. More...
 
SelfSmoothingSigmasAreSpecifiedInPhysicalUnitsOn ()
 Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels. More...
 
SelfSmoothingSigmasAreSpecifiedInPhysicalUnitsOff ()
 Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels. More...
 
- Public Member Functions inherited from itk::simple::ProcessObject
virtual void Abort ()
 
virtual int AddCommand (itk::simple::EventEnum event, const std::function< void()> &func)
 Directly add a callback to observe an event. More...
 
virtual int AddCommand (itk::simple::EventEnum event, itk::simple::Command &cmd)
 Add a Command Object to observer the event. More...
 
virtual float GetProgress () const
 An Active Measurement of the progress of execution. More...
 
virtual bool HasCommand (itk::simple::EventEnum event) const
 Query of this object has any registered commands for event. More...
 
 ProcessObject ()
 
virtual void RemoveAllCommands ()
 Remove all registered commands. More...
 
virtual ~ProcessObject ()
 
virtual void DebugOn ()
 
virtual void DebugOff ()
 
virtual bool GetDebug () const
 
virtual void SetDebug (bool debugFlag)
 
virtual void SetNumberOfThreads (unsigned int n)
 
virtual unsigned int GetNumberOfThreads () const
 
virtual void SetNumberOfWorkUnits (unsigned int n)
 
virtual unsigned int GetNumberOfWorkUnits () const
 

Protected Member Functions

unsigned long AddITKObserver (const itk::EventObject &, itk::Command *) override
 
template<class TImageType >
itk::ImageToImageMetricv4< TImageType, TImageType, TImageType, double, itk::DefaultImageToImageMetricTraitsv4< TImageType, TImageType, TImageType, double > > * CreateMetric ()
 
itk::ObjectToObjectOptimizerBaseTemplate< double > * CreateOptimizer (unsigned int numberOfTransformParameters)
 
template<typename TMetric >
itk::RegistrationParameterScalesEstimator< TMetric > * CreateScalesEstimator ()
 
template<unsigned int VDimension>
itk::SpatialObject< VDimension > * CreateSpatialObjectMask (const Image &mask)
 
template<typename TTransformAdaptorPointer , typename TRegistrationMethod >
std::vector< TTransformAdaptorPointer > CreateTransformParametersAdaptor (TRegistrationMethod *method)
 
template<class TImage >
double EvaluateInternal (const Image &fixed, const Image &moving)
 
template<class TImage >
Transform ExecuteInternal (const Image &fixed, const Image &moving)
 
void OnActiveProcessDelete () noexcept override
 
void PreUpdate (itk::ProcessObject *p) override
 
void RemoveITKObserver (EventCommand &e) override
 
template<class TImageType >
void SetupMetric (itk::ImageToImageMetricv4< TImageType, TImageType, TImageType, double, itk::DefaultImageToImageMetricTraitsv4< TImageType, TImageType, TImageType, double > > *, const TImageType *, const TImageType *)
 
- Protected Member Functions inherited from itk::simple::ProcessObject
virtual itk::ProcessObjectGetActiveProcess ()
 
virtual void onCommandDelete (const itk::simple::Command *cmd) noexcept
 
- Protected Member Functions inherited from itk::simple::NonCopyable
 NonCopyable ()=default
 
 NonCopyable (const NonCopyable &)=delete
 
NonCopyableoperator= (const NonCopyable &)=delete
 

Private Types

typedef double(ImageRegistrationMethod::* EvaluateMemberFunctionType) (const Image &fixed, const Image &moving)
 
typedef Transform(ImageRegistrationMethod::* MemberFunctionType) (const Image &fixed, const Image &moving)
 
enum  MetricType {
  ANTSNeighborhoodCorrelation,
  Correlation,
  Demons,
  JointHistogramMutualInformation,
  MeanSquares,
  MattesMutualInformation
}
 
enum  OptimizerScalesType {
  Manual,
  Jacobian,
  IndexShift,
  PhysicalShift
}
 
enum  OptimizerType {
  ConjugateGradientLineSearch,
  RegularStepGradientDescent,
  GradientDescent,
  GradientDescentLineSearch,
  LBFGSB,
  Exhaustive,
  Amoeba,
  Powell,
  OnePlusOneEvolutionary,
  LBFGS2
}
 

Private Attributes

itk::ObjectToObjectOptimizerBaseTemplate< double > * m_ActiveOptimizer
 
std::unique_ptr< detail::MemberFunctionFactory< EvaluateMemberFunctionType > > m_EvaluateMemberFactory
 
Transform m_FixedInitialTransform
 
Transform m_InitialTransform
 
bool m_InitialTransformInPlace
 
InterpolatorEnum m_Interpolator
 
unsigned int m_Iteration
 
std::unique_ptr< detail::MemberFunctionFactory< MemberFunctionType > > m_MemberFactory
 
Image m_MetricFixedMaskImage
 
double m_MetricIntensityDifferenceThreshold
 
Image m_MetricMovingMaskImage
 
unsigned int m_MetricNumberOfHistogramBins
 
unsigned int m_MetricRadius
 
std::vector< double > m_MetricSamplingPercentage
 
unsigned int m_MetricSamplingSeed
 
MetricSamplingStrategyType m_MetricSamplingStrategy
 
MetricType m_MetricType
 
bool m_MetricUseFixedImageGradientFilter
 
bool m_MetricUseMovingImageGradientFilter
 
double m_MetricValue
 
double m_MetricVarianceForJointPDFSmoothing
 
Transform m_MovingInitialTransform
 
uint64_t m_NumberOfValidPoints
 
double m_OptimizerConvergenceMinimumValue
 
unsigned int m_OptimizerConvergenceWindowSize
 
double m_OptimizerCostFunctionConvergenceFactor
 
unsigned int m_OptimizerDeltaConvergenceDistance
 
double m_OptimizerDeltaConvergenceTolerance
 
double m_OptimizerEpsilon
 
EstimateLearningRateType m_OptimizerEstimateLearningRate
 
double m_OptimizerFunctionConvergenceTolerance
 
double m_OptimizerGradientConvergenceTolerance
 
double m_OptimizerGradientMagnitudeTolerance
 
double m_OptimizerGrowthFactor
 
unsigned int m_OptimizerHessianApproximationAccuracy
 
double m_OptimizerInitialRadius
 
double m_OptimizerLearningRate
 
double m_OptimizerLineSearchAccuracy
 
double m_OptimizerLineSearchEpsilon
 
double m_OptimizerLineSearchLowerLimit
 
unsigned int m_OptimizerLineSearchMaximumEvaluations
 
unsigned int m_OptimizerLineSearchMaximumIterations
 
double m_OptimizerLineSearchMaximumStep
 
double m_OptimizerLineSearchMinimumStep
 
double m_OptimizerLineSearchUpperLimit
 
double m_OptimizerLowerBound
 
unsigned int m_OptimizerMaximumLineIterations
 
unsigned int m_OptimizerMaximumNumberOfCorrections
 
unsigned int m_OptimizerMaximumNumberOfFunctionEvaluations
 
double m_OptimizerMaximumStepSizeInPhysicalUnits
 
double m_OptimizerMinimumStepLength
 
unsigned int m_OptimizerNumberOfIterations
 
std::vector< unsigned int > m_OptimizerNumberOfSteps
 
double m_OptimizerParametersConvergenceTolerance
 
double m_OptimizerRelaxationFactor
 
std::vector< double > m_OptimizerScales
 
unsigned int m_OptimizerScalesCentralRegionRadius
 
double m_OptimizerScalesSmallParameterVariation
 
OptimizerScalesType m_OptimizerScalesType
 
unsigned int m_OptimizerSeed
 
double m_OptimizerShrinkFactor
 
double m_OptimizerSimplexDelta
 
double m_OptimizerSolutionAccuracy
 
double m_OptimizerStepLength
 
double m_OptimizerStepTolerance
 
bool m_OptimizerTrace
 
OptimizerType m_OptimizerType
 
double m_OptimizerUpperBound
 
double m_OptimizerValueTolerance
 
std::vector< double > m_OptimizerWeights
 
bool m_OptimizerWithRestarts
 
std::function< unsigned int()> m_pfGetCurrentLevel
 
std::function< uint64_t()> m_pfGetMetricNumberOfValidPoints
 
std::function< double()> m_pfGetMetricValue
 
std::function< double()> m_pfGetOptimizerConvergenceValue
 
std::function< unsigned int()> m_pfGetOptimizerIteration
 
std::function< double()> m_pfGetOptimizerLearningRate
 
std::function< std::vector< double >)> m_pfGetOptimizerPosition
 
std::function< std::vector< double >)> m_pfGetOptimizerScales
 
std::function< std::string()> m_pfGetOptimizerStopConditionDescription
 
std::function< void(itk::TransformBase *outTransform)> m_pfUpdateWithBestValue
 
std::vector< unsigned int > m_ShrinkFactorsPerLevel
 
bool m_SmoothingSigmasAreSpecifiedInPhysicalUnits
 
std::vector< double > m_SmoothingSigmasPerLevel
 
std::string m_StopConditionDescription
 
std::vector< unsigned int > m_TransformBSplineScaleFactors
 
std::vector< double > m_VirtualDomainDirection
 
std::vector< double > m_VirtualDomainOrigin
 
std::vector< uint32_tm_VirtualDomainSize
 
std::vector< double > m_VirtualDomainSpacing
 

Friends

struct detail::MemberFunctionAddressor< MemberFunctionType >
 

Additional Inherited Members

- Static Public Member Functions inherited from itk::simple::ProcessObject
static bool GetGlobalDefaultDebug ()
 
static void GlobalDefaultDebugOff ()
 
static void GlobalDefaultDebugOn ()
 
static void SetGlobalDefaultDebug (bool debugFlag)
 
static void GlobalWarningDisplayOn ()
 
static void GlobalWarningDisplayOff ()
 
static void SetGlobalWarningDisplay (bool flag)
 
static bool GetGlobalWarningDisplay ()
 
static double GetGlobalDefaultCoordinateTolerance ()
 Access the global tolerance to determine congruent spaces. More...
 
static void SetGlobalDefaultCoordinateTolerance (double)
 Access the global tolerance to determine congruent spaces. More...
 
static double GetGlobalDefaultDirectionTolerance ()
 Access the global tolerance to determine congruent spaces. More...
 
static void SetGlobalDefaultDirectionTolerance (double)
 Access the global tolerance to determine congruent spaces. More...
 
static bool SetGlobalDefaultThreader (const std::string &threader)
 Set/Get the default threader used for process objects. More...
 
static std::string GetGlobalDefaultThreader ()
 Set/Get the default threader used for process objects. More...
 
static void SetGlobalDefaultNumberOfThreads (unsigned int n)
 
static unsigned int GetGlobalDefaultNumberOfThreads ()
 Set/Get the default threader used for process objects. More...
 
- Static Protected Member Functions inherited from itk::simple::ProcessObject
template<class TImageType >
static TImageType::ConstPointer CastImageToITK (const Image &img)
 
template<class TPixelType , unsigned int VImageDimension, unsigned int VLength, template< typename, unsigned int > class TVector>
static Image CastITKToImage (itk::Image< TVector< TPixelType, VLength >, VImageDimension > *img)
 
template<class TImageType >
static Image CastITKToImage (TImageType *img)
 
static const itk::EventObjectGetITKEventObject (EventEnum e)
 
template<typename T >
static std::ostream & ToStringHelper (std::ostream &os, const T &v)
 
static std::ostream & ToStringHelper (std::ostream &os, const char &v)
 
static std::ostream & ToStringHelper (std::ostream &os, const signed char &v)
 
static std::ostream & ToStringHelper (std::ostream &os, const unsigned char &v)
 

Detailed Description

An interface method to the modular ITKv4 registration framework.

This interface method class encapsulates typical registration usage by incorporating all the necessary elements for performing a simple image registration between two images. This method also allows for multistage registration whereby each stage is characterized by possibly different transforms and different image metrics. For example, many users will want to perform a linear registration followed by deformable registration where both stages are performed in multiple levels. Each level can be characterized by:

Multiple stages are handled by linking multiple instantiations of this class where the output transform is added to the optional composite transform input.

See also
itk::ImageRegistrationMethodv4
itk::ImageToImageMetricv4
itk::ObjectToObjectOptimizerBaseTemplate
Examples
ImageRegistrationMethod1/ImageRegistrationMethod1.cxx, ImageRegistrationMethod2/ImageRegistrationMethod2.cxx, ImageRegistrationMethodBSpline1/ImageRegistrationMethodBSpline1.cxx, ImageRegistrationMethodBSpline3/ImageRegistrationMethodBSpline3.cxx, and ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx.

Definition at line 89 of file sitkImageRegistrationMethod.h.

Member Typedef Documentation

◆ EvaluateMemberFunctionType

typedef double(ImageRegistrationMethod::* itk::simple::ImageRegistrationMethod::EvaluateMemberFunctionType) (const Image &fixed, const Image &moving)
private

Definition at line 699 of file sitkImageRegistrationMethod.h.

◆ MemberFunctionType

typedef Transform(ImageRegistrationMethod::* itk::simple::ImageRegistrationMethod::MemberFunctionType) (const Image &fixed, const Image &moving)
private

Definition at line 698 of file sitkImageRegistrationMethod.h.

◆ Self

Definition at line 94 of file sitkImageRegistrationMethod.h.

◆ Superclass

Definition at line 95 of file sitkImageRegistrationMethod.h.

Member Enumeration Documentation

◆ EstimateLearningRateType

Enumerator
Never 

Never run estimation, use provided value.

Once 

Estimate learning once each level, ignore provided values.

EachIteration 

Estimate learning rate at each iteration.

Examples
ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx.

Definition at line 255 of file sitkImageRegistrationMethod.h.

◆ MetricSamplingStrategyType

Enumerator
NONE 
REGULAR 
RANDOM 

Definition at line 498 of file sitkImageRegistrationMethod.h.

◆ MetricType

Enumerator
ANTSNeighborhoodCorrelation 
Correlation 
Demons 
JointHistogramMutualInformation 
MeanSquares 
MattesMutualInformation 

Definition at line 788 of file sitkImageRegistrationMethod.h.

◆ OptimizerScalesType

Enumerator
Manual 
Jacobian 
IndexShift 
PhysicalShift 

Definition at line 776 of file sitkImageRegistrationMethod.h.

◆ OptimizerType

Enumerator
ConjugateGradientLineSearch 
RegularStepGradientDescent 
GradientDescent 
GradientDescentLineSearch 
LBFGSB 
Exhaustive 
Amoeba 
Powell 
OnePlusOneEvolutionary 
LBFGS2 

Definition at line 716 of file sitkImageRegistrationMethod.h.

Constructor & Destructor Documentation

◆ ~ImageRegistrationMethod()

itk::simple::ImageRegistrationMethod::~ImageRegistrationMethod ( )
override

◆ ImageRegistrationMethod()

itk::simple::ImageRegistrationMethod::ImageRegistrationMethod ( )

Member Function Documentation

◆ AddITKObserver()

unsigned long itk::simple::ImageRegistrationMethod::AddITKObserver ( const itk::EventObject ,
itk::Command  
)
overrideprotectedvirtual

Reimplemented from itk::simple::ProcessObject.

◆ CreateMetric()

template<class TImageType >
itk::ImageToImageMetricv4<TImageType, TImageType, TImageType, double, itk::DefaultImageToImageMetricTraitsv4< TImageType, TImageType, TImageType, double > >* itk::simple::ImageRegistrationMethod::CreateMetric ( )
protected

◆ CreateOptimizer()

itk::ObjectToObjectOptimizerBaseTemplate<double>* itk::simple::ImageRegistrationMethod::CreateOptimizer ( unsigned int  numberOfTransformParameters)
protected

◆ CreateScalesEstimator()

template<typename TMetric >
itk::RegistrationParameterScalesEstimator< TMetric >* itk::simple::ImageRegistrationMethod::CreateScalesEstimator ( )
protected

◆ CreateSpatialObjectMask()

template<unsigned int VDimension>
itk::SpatialObject<VDimension>* itk::simple::ImageRegistrationMethod::CreateSpatialObjectMask ( const Image mask)
protected

◆ CreateTransformParametersAdaptor()

template<typename TTransformAdaptorPointer , typename TRegistrationMethod >
std::vector< TTransformAdaptorPointer > itk::simple::ImageRegistrationMethod::CreateTransformParametersAdaptor ( TRegistrationMethod *  method)
protected

◆ EvaluateInternal()

template<class TImage >
double itk::simple::ImageRegistrationMethod::EvaluateInternal ( const Image fixed,
const Image moving 
)
protected

◆ Execute()

Transform itk::simple::ImageRegistrationMethod::Execute ( const Image fixed,
const Image moving 
)

◆ ExecuteInternal()

template<class TImage >
Transform itk::simple::ImageRegistrationMethod::ExecuteInternal ( const Image fixed,
const Image moving 
)
protected

◆ GetCurrentLevel()

unsigned int itk::simple::ImageRegistrationMethod::GetCurrentLevel ( ) const

◆ GetFixedInitialTransform()

Transform itk::simple::ImageRegistrationMethod::GetFixedInitialTransform ( ) const
inline

Set transform mapping to the fixed domain.

This transform is used to map from the virtual domain to the fixed image domain.

By default this transform is an identity, and not used.

See also
itk::ImageRegistrationMethodv4::SetFixedInitialTransform

Definition at line 198 of file sitkImageRegistrationMethod.h.

◆ GetInitialTransform()

Transform itk::simple::ImageRegistrationMethod::GetInitialTransform ( )
inline

Set the initial transform and parameters to optimize.

This transform is applied before the MovingInitialTransform, to map from the virtual image domain to the moving image domain.

If the inPlace flag is explicitly false, then the ITK registration will internally make a copy, and the transform will not be accessible during registration. Otherwise, the accessible InitialTransform value will be the same object used during registration, and will have a modified value upon completion.

See also
itk::ImageRegistrationMethodv4::SetInitialTransform

Definition at line 142 of file sitkImageRegistrationMethod.h.

◆ GetInitialTransformInPlace()

bool itk::simple::ImageRegistrationMethod::GetInitialTransformInPlace ( ) const
inline

Set the initial transform and parameters to optimize.

This transform is applied before the MovingInitialTransform, to map from the virtual image domain to the moving image domain.

If the inPlace flag is explicitly false, then the ITK registration will internally make a copy, and the transform will not be accessible during registration. Otherwise, the accessible InitialTransform value will be the same object used during registration, and will have a modified value upon completion.

See also
itk::ImageRegistrationMethodv4::SetInitialTransform

Definition at line 144 of file sitkImageRegistrationMethod.h.

◆ GetInterpolator()

InterpolatorEnum itk::simple::ImageRegistrationMethod::GetInterpolator ( )
inline

Set and get the interpolator to use.

Definition at line 115 of file sitkImageRegistrationMethod.h.

◆ GetMetricNumberOfValidPoints()

uint64_t itk::simple::ImageRegistrationMethod::GetMetricNumberOfValidPoints ( ) const

Current number of points used of metric evaluation

This is a active measurement connected to the registration processes during registration. This number is number of point in the virtual domain which overlap the fixed image and the moving image. It is valid for sparse or dense sampling. After execution of registration this will contain the last value.

◆ GetMetricSamplingPercentagePerLevel()

const std::vector<double>& itk::simple::ImageRegistrationMethod::GetMetricSamplingPercentagePerLevel ( ) const

Get the percentage of pixels used for metric evaluation.

◆ GetMetricValue()

double itk::simple::ImageRegistrationMethod::GetMetricValue ( ) const

◆ GetMovingInitialTransform()

Transform itk::simple::ImageRegistrationMethod::GetMovingInitialTransform ( ) const
inline

Set a fixed transform component towards moving domain.

The InitialTransform is added to this transform to form the composite transform which maps from the virtual domain to the moving image's domain. The parameters of this transform are not optimized. This is an advanced option used when the output transformed with be composed into a composite transform.

By default this transform is an identity, and not used.

See also
itk::ImageRegistrationMethodv4::SetMovingInitialTransform

Definition at line 182 of file sitkImageRegistrationMethod.h.

◆ GetName()

std::string itk::simple::ImageRegistrationMethod::GetName ( ) const
inlineoverridevirtual

return user readable name for the filter

Implements itk::simple::ProcessObject.

Definition at line 101 of file sitkImageRegistrationMethod.h.

◆ GetOptimizerConvergenceValue()

double itk::simple::ImageRegistrationMethod::GetOptimizerConvergenceValue ( ) const

◆ GetOptimizerIteration()

unsigned int itk::simple::ImageRegistrationMethod::GetOptimizerIteration ( ) const

◆ GetOptimizerLearningRate()

double itk::simple::ImageRegistrationMethod::GetOptimizerLearningRate ( ) const

◆ GetOptimizerPosition()

std::vector<double> itk::simple::ImageRegistrationMethod::GetOptimizerPosition ( ) const

◆ GetOptimizerScales()

std::vector<double> itk::simple::ImageRegistrationMethod::GetOptimizerScales ( ) const

Get the OptimizerScales.

If the scales are explicitly set then this method returns those values. If an estimator is used then this is an active measurement returning the scales estimated by the estimator and is only available during execution.

◆ GetOptimizerStopConditionDescription()

std::string itk::simple::ImageRegistrationMethod::GetOptimizerStopConditionDescription ( ) const

◆ GetOptimizerWeights()

std::vector<double> itk::simple::ImageRegistrationMethod::GetOptimizerWeights ( ) const

A per parameter weighting array for the optimizer.

Allows setting of a per-local-parameter weighting array. If unset, the weights are treated as identity. Weights are multiplied by the gradient at the same time scaling is applied. Weights are similar to the scales but not estimated, and may be used, for example, to easily mask out a particular parameter during optimization to hold it constant. Or they may be used to apply another kind of prior knowledge.

◆ MetricEvaluate()

double itk::simple::ImageRegistrationMethod::MetricEvaluate ( const Image fixed,
const Image moving 
)

Get the value of the metric given the state of the method.

Passing a fixed and moving image, this method constructs and configures a metric object to obtain the value. This will take into consideration the current transforms, metric, interpolator, and image masks. It does not take into consideration the sampling strategy, smoothing sigmas, or the shrink factors.

◆ MetricUseFixedImageGradientFilterOff()

Self& itk::simple::ImageRegistrationMethod::MetricUseFixedImageGradientFilterOff ( )
inline

Enable image gradient computation by a filter.

By default the image gradient is computed by itk::GradientRecursiveGaussianImageFiter. If disabled then a central difference function with be computed as needed.

See also
itk::ImageToImageMetricv4::SetUseFixedImageGradientFilter
Examples
ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx.

Definition at line 521 of file sitkImageRegistrationMethod.h.

◆ MetricUseFixedImageGradientFilterOn()

Self& itk::simple::ImageRegistrationMethod::MetricUseFixedImageGradientFilterOn ( )
inline

Enable image gradient computation by a filter.

By default the image gradient is computed by itk::GradientRecursiveGaussianImageFiter. If disabled then a central difference function with be computed as needed.

See also
itk::ImageToImageMetricv4::SetUseFixedImageGradientFilter

Definition at line 520 of file sitkImageRegistrationMethod.h.

◆ MetricUseMovingImageGradientFilterOff()

Self& itk::simple::ImageRegistrationMethod::MetricUseMovingImageGradientFilterOff ( )
inline

Enable image gradient computation by a filter.

By default the image gradient is computed by itk::GradientRecursiveGaussianImageFiter. If disabled then a central difference function will be computed for each sample as needed.

See also
itk::ImageToImageMetricv4::SetUseMovingImageGradientFilter

Definition at line 537 of file sitkImageRegistrationMethod.h.

◆ MetricUseMovingImageGradientFilterOn()

Self& itk::simple::ImageRegistrationMethod::MetricUseMovingImageGradientFilterOn ( )
inline

Enable image gradient computation by a filter.

By default the image gradient is computed by itk::GradientRecursiveGaussianImageFiter. If disabled then a central difference function will be computed for each sample as needed.

See also
itk::ImageToImageMetricv4::SetUseMovingImageGradientFilter

Definition at line 536 of file sitkImageRegistrationMethod.h.

◆ OnActiveProcessDelete()

void itk::simple::ImageRegistrationMethod::OnActiveProcessDelete ( )
overrideprotectedvirtualnoexcept

Reimplemented from itk::simple::ProcessObject.

◆ PreUpdate()

void itk::simple::ImageRegistrationMethod::PreUpdate ( itk::ProcessObject p)
overrideprotectedvirtual

Reimplemented from itk::simple::ProcessObject.

◆ RemoveITKObserver()

void itk::simple::ImageRegistrationMethod::RemoveITKObserver ( EventCommand e)
overrideprotectedvirtual

Reimplemented from itk::simple::ProcessObject.

◆ SetFixedInitialTransform()

Self& itk::simple::ImageRegistrationMethod::SetFixedInitialTransform ( const Transform transform)
inline

Set transform mapping to the fixed domain.

This transform is used to map from the virtual domain to the fixed image domain.

By default this transform is an identity, and not used.

See also
itk::ImageRegistrationMethodv4::SetFixedInitialTransform

Definition at line 196 of file sitkImageRegistrationMethod.h.

◆ SetInitialTransform() [1/2]

Self& itk::simple::ImageRegistrationMethod::SetInitialTransform ( const Transform transform)

Set the initial transform and parameters to optimize.

This transform is applied before the MovingInitialTransform, to map from the virtual image domain to the moving image domain.

If the inPlace flag is explicitly false, then the ITK registration will internally make a copy, and the transform will not be accessible during registration. Otherwise, the accessible InitialTransform value will be the same object used during registration, and will have a modified value upon completion.

See also
itk::ImageRegistrationMethodv4::SetInitialTransform
Examples
ImageRegistrationMethod1/ImageRegistrationMethod1.cxx, ImageRegistrationMethod2/ImageRegistrationMethod2.cxx, ImageRegistrationMethodBSpline1/ImageRegistrationMethodBSpline1.cxx, and ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx.

◆ SetInitialTransform() [2/2]

Self& itk::simple::ImageRegistrationMethod::SetInitialTransform ( Transform transform,
bool  inPlace = true 
)

Set the initial transform and parameters to optimize.

This transform is applied before the MovingInitialTransform, to map from the virtual image domain to the moving image domain.

If the inPlace flag is explicitly false, then the ITK registration will internally make a copy, and the transform will not be accessible during registration. Otherwise, the accessible InitialTransform value will be the same object used during registration, and will have a modified value upon completion.

See also
itk::ImageRegistrationMethodv4::SetInitialTransform

◆ SetInitialTransformAsBSpline()

void itk::simple::ImageRegistrationMethod::SetInitialTransformAsBSpline ( BSplineTransform transform,
bool  inPlace = true,
const std::vector< unsigned int > &  scaleFactors = std::vector< unsigned int >() 
)

Set an initial BSpline transform to optimize.

A specialization of SetInitialTransform for BSplineTransforms which can take an additional scaleFactors parameter. The scaleFactors specifies the a isotropic scaling factor per level for the BSpline transform mesh size with respect to the initial transform. For example to double the BSpline mesh resolution at each of 3 levels the vector [1,2,4] should be provided.

If a per level scale factor is 0 or omitted than no transform adapter will be created for that level.

See also
itk::BSplineTransformParametersAdaptor
Examples
ImageRegistrationMethodBSpline3/ImageRegistrationMethodBSpline3.cxx.

◆ SetInterpolator()

Self& itk::simple::ImageRegistrationMethod::SetInterpolator ( InterpolatorEnum  Interpolator)
inline

◆ SetMetricAsANTSNeighborhoodCorrelation()

Self& itk::simple::ImageRegistrationMethod::SetMetricAsANTSNeighborhoodCorrelation ( unsigned int  radius)

Use normalized cross correlation using a small neighborhood for each voxel between two images, with speed optimizations for dense registration.

See also
itk::ANTSNeighborhoodCorrelationImageToImageMetricv4
Examples
ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx.

◆ SetMetricAsCorrelation()

Self& itk::simple::ImageRegistrationMethod::SetMetricAsCorrelation ( )

Use negative normalized cross correlation image metric.

See also
itk::CorrelationImageToImageMetricv4
Examples
ImageRegistrationMethodBSpline1/ImageRegistrationMethodBSpline1.cxx.

◆ SetMetricAsDemons()

Self& itk::simple::ImageRegistrationMethod::SetMetricAsDemons ( double  intensityDifferenceThreshold = 0.001)

Use demons image metric.

See also
itk::DemonsImageToImageMetricv4

◆ SetMetricAsJointHistogramMutualInformation()

Self& itk::simple::ImageRegistrationMethod::SetMetricAsJointHistogramMutualInformation ( unsigned int  numberOfHistogramBins = 20,
double  varianceForJointPDFSmoothing = 1.5 
)

◆ SetMetricAsMattesMutualInformation()

Self& itk::simple::ImageRegistrationMethod::SetMetricAsMattesMutualInformation ( unsigned int  numberOfHistogramBins = 50)

Use the mutual information between two images to be registered using the method of Mattes et al.

See also
itk::MattesMutualInformationImageToImageMetricv4

◆ SetMetricAsMeanSquares()

Self& itk::simple::ImageRegistrationMethod::SetMetricAsMeanSquares ( )

Use negative means squares image metric.

See also
itk::MeanSquaresImageToImageMetricv4
Examples
ImageRegistrationMethod1/ImageRegistrationMethod1.cxx.

◆ SetMetricFixedMask()

Self& itk::simple::ImageRegistrationMethod::SetMetricFixedMask ( const Image binaryMask)

Set an image mask in order to restrict the sampled points for the metric.

The image is expected to be in the same physical space as the FixedImage, and if the pixel type is not UInt8 than the image will base cast.

See also
itk::ImageToImageMetricv4::SetFixedImageMask

◆ SetMetricMovingMask()

Self& itk::simple::ImageRegistrationMethod::SetMetricMovingMask ( const Image binaryMask)

Set an image mask in order to restrict the sampled points for the metric in the moving image space.

The image is expected to be in the same physical space as the MovingImage, and if the pixel type is not UInt8 than the image will base cast.

See also
itk::ImageToImageMetricv4::SetMovingImageMask

◆ SetMetricSamplingPercentage()

Self& itk::simple::ImageRegistrationMethod::SetMetricSamplingPercentage ( double  percentage,
unsigned int  seed = sitkWallClock 
)

Set percentage of pixels sampled for metric evaluation.

The percentage is of the number of pixels in the virtual domain per level after the shrink factor has been applied.

The seed parameter is used to seed the pseudo-random number generator used in generating the sampling set of points. If the seed parameter is 0, then the wall clock is used, otherwise the fixed seed is used for reproducible behavior.

When providing multiple values the number of entries must match the number of shrink factors and smoothing sigmas. If a single value is given, then the same percentage is used for all levels.

See also
itk::ImageRegistrationMethodv4::SetMetricSamplingPercentage

◆ SetMetricSamplingPercentagePerLevel()

Self& itk::simple::ImageRegistrationMethod::SetMetricSamplingPercentagePerLevel ( const std::vector< double > &  percentage,
unsigned int  seed = sitkWallClock 
)

Set percentage of pixels sampled for metric evaluation.

The percentage is of the number of pixels in the virtual domain per level after the shrink factor has been applied.

The seed parameter is used to seed the pseudo-random number generator used in generating the sampling set of points. If the seed parameter is 0, then the wall clock is used, otherwise the fixed seed is used for reproducible behavior.

When providing multiple values the number of entries must match the number of shrink factors and smoothing sigmas. If a single value is given, then the same percentage is used for all levels.

See also
itk::ImageRegistrationMethodv4::SetMetricSamplingPercentage

◆ SetMetricSamplingStrategy()

Self& itk::simple::ImageRegistrationMethod::SetMetricSamplingStrategy ( MetricSamplingStrategyType  strategy)

Set sampling strategy for sample generation.

See also
itk::ImageRegistrationMethodv4::SetMetricSamplingStrategy

◆ SetMetricUseFixedImageGradientFilter()

Self& itk::simple::ImageRegistrationMethod::SetMetricUseFixedImageGradientFilter ( bool  )

Enable image gradient computation by a filter.

By default the image gradient is computed by itk::GradientRecursiveGaussianImageFiter. If disabled then a central difference function with be computed as needed.

See also
itk::ImageToImageMetricv4::SetUseFixedImageGradientFilter

◆ SetMetricUseMovingImageGradientFilter()

Self& itk::simple::ImageRegistrationMethod::SetMetricUseMovingImageGradientFilter ( bool  )

Enable image gradient computation by a filter.

By default the image gradient is computed by itk::GradientRecursiveGaussianImageFiter. If disabled then a central difference function will be computed for each sample as needed.

See also
itk::ImageToImageMetricv4::SetUseMovingImageGradientFilter

◆ SetMovingInitialTransform()

Self& itk::simple::ImageRegistrationMethod::SetMovingInitialTransform ( const Transform transform)
inline

Set a fixed transform component towards moving domain.

The InitialTransform is added to this transform to form the composite transform which maps from the virtual domain to the moving image's domain. The parameters of this transform are not optimized. This is an advanced option used when the output transformed with be composed into a composite transform.

By default this transform is an identity, and not used.

See also
itk::ImageRegistrationMethodv4::SetMovingInitialTransform
Examples
ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx.

Definition at line 180 of file sitkImageRegistrationMethod.h.

◆ SetOptimizerAsAmoeba()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerAsAmoeba ( double  simplexDelta,
unsigned int  numberOfIterations,
double  parametersConvergenceTolerance = 1e-8,
double  functionConvergenceTolerance = 1e-4,
bool  withRestarts = false 
)

Set optimizer to Nelder-Mead downhill simplex algorithm.

See also
itk::AmoebaOptimizerv4

◆ SetOptimizerAsConjugateGradientLineSearch()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerAsConjugateGradientLineSearch ( double  learningRate,
unsigned int  numberOfIterations,
double  convergenceMinimumValue = 1e-6,
unsigned int  convergenceWindowSize = 10,
double  lineSearchLowerLimit = 0,
double  lineSearchUpperLimit = 5.0,
double  lineSearchEpsilon = 0.01,
unsigned int  lineSearchMaximumIterations = 20,
EstimateLearningRateType  estimateLearningRate = Once,
double  maximumStepSizeInPhysicalUnits = 0.0 
)

Conjugate gradient descent optimizer with a golden section line search for nonlinear optimization.

See also
itk::ConjugateGradientLineSearchOptimizerv4Template

◆ SetOptimizerAsExhaustive()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerAsExhaustive ( const std::vector< unsigned int > &  numberOfSteps,
double  stepLength = 1.0 
)

Set the optimizer to sample the metric at regular steps.

At each iteration the GetOptimizerIteration, can be used to index into the sampling grid along with the GetCurrentMetricValue.

The resulting transform and value at the end of execution is the best location.

The OptimizerScales can be used to perform anisotropic sampling.

Note
This optimizer is not suitable for use in conjunction with the multiple scales.
See also
itk::ExhaustiveOptimizerv4

◆ SetOptimizerAsGradientDescent()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerAsGradientDescent ( double  learningRate,
unsigned int  numberOfIterations,
double  convergenceMinimumValue = 1e-6,
unsigned int  convergenceWindowSize = 10,
EstimateLearningRateType  estimateLearningRate = Once,
double  maximumStepSizeInPhysicalUnits = 0.0 
)

◆ SetOptimizerAsGradientDescentLineSearch()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerAsGradientDescentLineSearch ( double  learningRate,
unsigned int  numberOfIterations,
double  convergenceMinimumValue = 1e-6,
unsigned int  convergenceWindowSize = 10,
double  lineSearchLowerLimit = 0,
double  lineSearchUpperLimit = 5.0,
double  lineSearchEpsilon = 0.01,
unsigned int  lineSearchMaximumIterations = 20,
EstimateLearningRateType  estimateLearningRate = Once,
double  maximumStepSizeInPhysicalUnits = 0.0 
)

◆ SetOptimizerAsLBFGS2()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerAsLBFGS2 ( double  solutionAccuracy = 1e-5,
unsigned int  numberOfIterations = 0,
unsigned int  hessianApproximateAccuracy = 6,
unsigned int  deltaConvergenceDistance = 0,
double  deltaConvergenceTolerance = 1e-5,
unsigned int  lineSearchMaximumEvaluations = 40,
double  lineSearchMinimumStep = 1e-20,
double  lineSearchMaximumStep = 1e20,
double  lineSearchAccuracy = 1e-4 
)

Limited memory Broyden Fletcher Goldfarb Shannon minimization without bounds.

The default parameters utilize LBFGSB in unbounded mode. This version is from LibLBFGS.

There are upto 3 stopping criteria:

  • the solution accuracy which is the magnitude of the gradient
  • the delta convergence which ensures the decrease of the metric
  • maximum number of iterations
See also
itk::LBFGS2Optimizerv4

◆ SetOptimizerAsLBFGSB()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerAsLBFGSB ( double  gradientConvergenceTolerance = 1e-5,
unsigned int  numberOfIterations = 500,
unsigned int  maximumNumberOfCorrections = 5,
unsigned int  maximumNumberOfFunctionEvaluations = 2000,
double  costFunctionConvergenceFactor = 1e+7,
double  lowerBound = std::numeric_limits< double >::min(),
double  upperBound = std::numeric_limits< double >::max(),
bool  trace = false 
)

Limited memory Broyden Fletcher Goldfarb Shannon minimization with simple bounds.

The default parameters utilize LBFGSB in unbounded mode.

See also
itk::LBFGSBOptimizerv4
Examples
ImageRegistrationMethodBSpline1/ImageRegistrationMethodBSpline1.cxx.

◆ SetOptimizerAsOnePlusOneEvolutionary()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerAsOnePlusOneEvolutionary ( unsigned int  numberOfIterations = 100,
double  epsilon = 1.5e-4,
double  initialRadius = 1.01,
double  growthFactor = -1.0,
double  shrinkFactor = -1.0,
unsigned int  seed = sitkWallClock 
)

1+1 evolutionary optimizer strategy.

The seed parameter is used to seed the pseudo-random number generator. If the seed parameter is 0, then the wall clock is used to seed, otherwise the fixed seed is used for reproducible behavior.

See also
itk::OnePlusOneEvolutionaryOptimizerv4

◆ SetOptimizerAsPowell()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerAsPowell ( unsigned int  numberOfIterations = 100,
unsigned int  maximumLineIterations = 100,
double  stepLength = 1,
double  stepTolerance = 1e-6,
double  valueTolerance = 1e-6 
)

Powell optimization using Brent line search.

See also
itk::PowellOptimizerv4

◆ SetOptimizerAsRegularStepGradientDescent()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerAsRegularStepGradientDescent ( double  learningRate,
double  minStep,
unsigned int  numberOfIterations,
double  relaxationFactor = 0.5,
double  gradientMagnitudeTolerance = 1e-4,
EstimateLearningRateType  estimateLearningRate = Never,
double  maximumStepSizeInPhysicalUnits = 0.0 
)

◆ SetOptimizerScales()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerScales ( const std::vector< double > &  scales)

Manually set per parameter weighting for the transform parameters.

◆ SetOptimizerScalesFromIndexShift()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerScalesFromIndexShift ( unsigned int  centralRegionRadius = 5,
double  smallParameterVariation = 0.01 
)

Estimate scales from maximum voxel shift in index space cause by parameter change.

See also
itk::RegistrationParameterScalesFromIndexShift

◆ SetOptimizerScalesFromJacobian()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerScalesFromJacobian ( unsigned int  centralRegionRadius = 5)

Estimate scales from Jacobian norms.

This scales estimator works well with versor based transforms.

See also
itk::RegistrationParameterScalesFromJacobian

◆ SetOptimizerScalesFromPhysicalShift()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerScalesFromPhysicalShift ( unsigned int  centralRegionRadius = 5,
double  smallParameterVariation = 0.01 
)

Estimating scales of transform parameters a step sizes, from the maximum voxel shift in physical space caused by a parameter change.

See also
itk::RegistrationParameterScalesFromPhysicalShift
Examples
ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx.

◆ SetOptimizerWeights()

Self& itk::simple::ImageRegistrationMethod::SetOptimizerWeights ( const std::vector< double > &  weights)

A per parameter weighting array for the optimizer.

Allows setting of a per-local-parameter weighting array. If unset, the weights are treated as identity. Weights are multiplied by the gradient at the same time scaling is applied. Weights are similar to the scales but not estimated, and may be used, for example, to easily mask out a particular parameter during optimization to hold it constant. Or they may be used to apply another kind of prior knowledge.

◆ SetShrinkFactorsPerLevel()

Self& itk::simple::ImageRegistrationMethod::SetShrinkFactorsPerLevel ( const std::vector< unsigned int > &  shrinkFactors)

Set the isotropic shrink factors for each level.

The virtual domain image is shrunk by this factor relative to the full size of the original virtual domain.

See also
itk::ImageRegistrationMethodv4::SetShrinkFactorsPerLevel
Examples
ImageRegistrationMethodBSpline3/ImageRegistrationMethodBSpline3.cxx, and ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx.

◆ SetSmoothingSigmasAreSpecifiedInPhysicalUnits()

Self& itk::simple::ImageRegistrationMethod::SetSmoothingSigmasAreSpecifiedInPhysicalUnits ( bool  arg)

Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels.

See also
itk::ImageRegistrationMethodv4::SetSmoothingSigmasAreSpecifiedInPhysicalUnits

◆ SetSmoothingSigmasPerLevel()

Self& itk::simple::ImageRegistrationMethod::SetSmoothingSigmasPerLevel ( const std::vector< double > &  smoothingSigmas)

Set the sigmas of Gaussian used for smoothing.

The smoothing is applied to both the fixed and the moving images at each level. The number of smoothing sigmas must match the number of shrink factors.

See also
itk::ImageRegistrationMethodv4::SetSmoothingSigmasPerLevel
Examples
ImageRegistrationMethodBSpline3/ImageRegistrationMethodBSpline3.cxx, and ImageRegistrationMethodDisplacement1/ImageRegistrationMethodDisplacement1.cxx.

◆ SetupMetric()

template<class TImageType >
void itk::simple::ImageRegistrationMethod::SetupMetric ( itk::ImageToImageMetricv4< TImageType, TImageType, TImageType, double, itk::DefaultImageToImageMetricTraitsv4< TImageType, TImageType, TImageType, double > > *  ,
const TImageType *  ,
const TImageType *   
)
protected

◆ SetVirtualDomain()

Self& itk::simple::ImageRegistrationMethod::SetVirtualDomain ( const std::vector< uint32_t > &  virtualSize,
const std::vector< double > &  virtualOrigin,
const std::vector< double > &  virtualSpacing,
const std::vector< double > &  virtualDirection 
)

Set the virtual domain used for sampling.

◆ SetVirtualDomainFromImage()

Self& itk::simple::ImageRegistrationMethod::SetVirtualDomainFromImage ( const Image virtualImage)

Set the virtual domain used for sampling.

◆ SmoothingSigmasAreSpecifiedInPhysicalUnitsOff()

Self& itk::simple::ImageRegistrationMethod::SmoothingSigmasAreSpecifiedInPhysicalUnitsOff ( )
inline

Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels.

See also
itk::ImageRegistrationMethodv4::SetSmoothingSigmasAreSpecifiedInPhysicalUnits

Definition at line 568 of file sitkImageRegistrationMethod.h.

◆ SmoothingSigmasAreSpecifiedInPhysicalUnitsOn()

Self& itk::simple::ImageRegistrationMethod::SmoothingSigmasAreSpecifiedInPhysicalUnitsOn ( )
inline

Enable the smoothing sigmas for each level in physical units (default) or in terms of voxels.

See also
itk::ImageRegistrationMethodv4::SetSmoothingSigmasAreSpecifiedInPhysicalUnits

Definition at line 567 of file sitkImageRegistrationMethod.h.

◆ ToString()

std::string itk::simple::ImageRegistrationMethod::ToString ( ) const
overridevirtual

Print the information about the object to a string.

If called when the process is being executed ( during a callback ), the ITK Optimizer and Transform objects will be printed.

Reimplemented from itk::simple::ProcessObject.

Friends And Related Function Documentation

◆ detail::MemberFunctionAddressor< MemberFunctionType >

Definition at line 700 of file sitkImageRegistrationMethod.h.

Member Data Documentation

◆ m_ActiveOptimizer

itk::ObjectToObjectOptimizerBaseTemplate<double>* itk::simple::ImageRegistrationMethod::m_ActiveOptimizer
private

Definition at line 820 of file sitkImageRegistrationMethod.h.

◆ m_EvaluateMemberFactory

std::unique_ptr<detail::MemberFunctionFactory<EvaluateMemberFunctionType> > itk::simple::ImageRegistrationMethod::m_EvaluateMemberFactory
private

Definition at line 702 of file sitkImageRegistrationMethod.h.

◆ m_FixedInitialTransform

Transform itk::simple::ImageRegistrationMethod::m_FixedInitialTransform
private

Definition at line 708 of file sitkImageRegistrationMethod.h.

◆ m_InitialTransform

Transform itk::simple::ImageRegistrationMethod::m_InitialTransform
private

Definition at line 705 of file sitkImageRegistrationMethod.h.

◆ m_InitialTransformInPlace

bool itk::simple::ImageRegistrationMethod::m_InitialTransformInPlace
private

Definition at line 706 of file sitkImageRegistrationMethod.h.

◆ m_Interpolator

InterpolatorEnum itk::simple::ImageRegistrationMethod::m_Interpolator
private

Definition at line 704 of file sitkImageRegistrationMethod.h.

◆ m_Iteration

unsigned int itk::simple::ImageRegistrationMethod::m_Iteration
private

Definition at line 817 of file sitkImageRegistrationMethod.h.

◆ m_MemberFactory

std::unique_ptr<detail::MemberFunctionFactory<MemberFunctionType> > itk::simple::ImageRegistrationMethod::m_MemberFactory
private

Definition at line 701 of file sitkImageRegistrationMethod.h.

◆ m_MetricFixedMaskImage

Image itk::simple::ImageRegistrationMethod::m_MetricFixedMaskImage
private

Definition at line 801 of file sitkImageRegistrationMethod.h.

◆ m_MetricIntensityDifferenceThreshold

double itk::simple::ImageRegistrationMethod::m_MetricIntensityDifferenceThreshold
private

Definition at line 797 of file sitkImageRegistrationMethod.h.

◆ m_MetricMovingMaskImage

Image itk::simple::ImageRegistrationMethod::m_MetricMovingMaskImage
private

Definition at line 802 of file sitkImageRegistrationMethod.h.

◆ m_MetricNumberOfHistogramBins

unsigned int itk::simple::ImageRegistrationMethod::m_MetricNumberOfHistogramBins
private

Definition at line 798 of file sitkImageRegistrationMethod.h.

◆ m_MetricRadius

unsigned int itk::simple::ImageRegistrationMethod::m_MetricRadius
private

Definition at line 796 of file sitkImageRegistrationMethod.h.

◆ m_MetricSamplingPercentage

std::vector<double> itk::simple::ImageRegistrationMethod::m_MetricSamplingPercentage
private

Definition at line 804 of file sitkImageRegistrationMethod.h.

◆ m_MetricSamplingSeed

unsigned int itk::simple::ImageRegistrationMethod::m_MetricSamplingSeed
private

Definition at line 806 of file sitkImageRegistrationMethod.h.

◆ m_MetricSamplingStrategy

MetricSamplingStrategyType itk::simple::ImageRegistrationMethod::m_MetricSamplingStrategy
private

Definition at line 805 of file sitkImageRegistrationMethod.h.

◆ m_MetricType

MetricType itk::simple::ImageRegistrationMethod::m_MetricType
private

Definition at line 795 of file sitkImageRegistrationMethod.h.

◆ m_MetricUseFixedImageGradientFilter

bool itk::simple::ImageRegistrationMethod::m_MetricUseFixedImageGradientFilter
private

Definition at line 808 of file sitkImageRegistrationMethod.h.

◆ m_MetricUseMovingImageGradientFilter

bool itk::simple::ImageRegistrationMethod::m_MetricUseMovingImageGradientFilter
private

Definition at line 809 of file sitkImageRegistrationMethod.h.

◆ m_MetricValue

double itk::simple::ImageRegistrationMethod::m_MetricValue
private

Definition at line 816 of file sitkImageRegistrationMethod.h.

◆ m_MetricVarianceForJointPDFSmoothing

double itk::simple::ImageRegistrationMethod::m_MetricVarianceForJointPDFSmoothing
private

Definition at line 799 of file sitkImageRegistrationMethod.h.

◆ m_MovingInitialTransform

Transform itk::simple::ImageRegistrationMethod::m_MovingInitialTransform
private

Definition at line 707 of file sitkImageRegistrationMethod.h.

◆ m_NumberOfValidPoints

uint64_t itk::simple::ImageRegistrationMethod::m_NumberOfValidPoints
private

Definition at line 818 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerConvergenceMinimumValue

double itk::simple::ImageRegistrationMethod::m_OptimizerConvergenceMinimumValue
private

Definition at line 739 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerConvergenceWindowSize

unsigned int itk::simple::ImageRegistrationMethod::m_OptimizerConvergenceWindowSize
private

Definition at line 740 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerCostFunctionConvergenceFactor

double itk::simple::ImageRegistrationMethod::m_OptimizerCostFunctionConvergenceFactor
private

Definition at line 744 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerDeltaConvergenceDistance

unsigned int itk::simple::ImageRegistrationMethod::m_OptimizerDeltaConvergenceDistance
private

Definition at line 764 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerDeltaConvergenceTolerance

double itk::simple::ImageRegistrationMethod::m_OptimizerDeltaConvergenceTolerance
private

Definition at line 765 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerEpsilon

double itk::simple::ImageRegistrationMethod::m_OptimizerEpsilon
private

Definition at line 757 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerEstimateLearningRate

EstimateLearningRateType itk::simple::ImageRegistrationMethod::m_OptimizerEstimateLearningRate
private

Definition at line 735 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerFunctionConvergenceTolerance

double itk::simple::ImageRegistrationMethod::m_OptimizerFunctionConvergenceTolerance
private

Definition at line 752 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerGradientConvergenceTolerance

double itk::simple::ImageRegistrationMethod::m_OptimizerGradientConvergenceTolerance
private

Definition at line 741 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerGradientMagnitudeTolerance

double itk::simple::ImageRegistrationMethod::m_OptimizerGradientMagnitudeTolerance
private

Definition at line 738 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerGrowthFactor

double itk::simple::ImageRegistrationMethod::m_OptimizerGrowthFactor
private

Definition at line 759 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerHessianApproximationAccuracy

unsigned int itk::simple::ImageRegistrationMethod::m_OptimizerHessianApproximationAccuracy
private

Definition at line 763 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerInitialRadius

double itk::simple::ImageRegistrationMethod::m_OptimizerInitialRadius
private

Definition at line 758 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerLearningRate

double itk::simple::ImageRegistrationMethod::m_OptimizerLearningRate
private

Definition at line 728 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerLineSearchAccuracy

double itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchAccuracy
private

Definition at line 769 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerLineSearchEpsilon

double itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchEpsilon
private

Definition at line 733 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerLineSearchLowerLimit

double itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchLowerLimit
private

Definition at line 731 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerLineSearchMaximumEvaluations

unsigned int itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchMaximumEvaluations
private

Definition at line 766 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerLineSearchMaximumIterations

unsigned int itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchMaximumIterations
private

Definition at line 734 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerLineSearchMaximumStep

double itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchMaximumStep
private

Definition at line 768 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerLineSearchMinimumStep

double itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchMinimumStep
private

Definition at line 767 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerLineSearchUpperLimit

double itk::simple::ImageRegistrationMethod::m_OptimizerLineSearchUpperLimit
private

Definition at line 732 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerLowerBound

double itk::simple::ImageRegistrationMethod::m_OptimizerLowerBound
private

Definition at line 745 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerMaximumLineIterations

unsigned int itk::simple::ImageRegistrationMethod::m_OptimizerMaximumLineIterations
private

Definition at line 754 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerMaximumNumberOfCorrections

unsigned int itk::simple::ImageRegistrationMethod::m_OptimizerMaximumNumberOfCorrections
private

Definition at line 742 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerMaximumNumberOfFunctionEvaluations

unsigned int itk::simple::ImageRegistrationMethod::m_OptimizerMaximumNumberOfFunctionEvaluations
private

Definition at line 743 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerMaximumStepSizeInPhysicalUnits

double itk::simple::ImageRegistrationMethod::m_OptimizerMaximumStepSizeInPhysicalUnits
private

Definition at line 736 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerMinimumStepLength

double itk::simple::ImageRegistrationMethod::m_OptimizerMinimumStepLength
private

Definition at line 729 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerNumberOfIterations

unsigned int itk::simple::ImageRegistrationMethod::m_OptimizerNumberOfIterations
private

Definition at line 730 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerNumberOfSteps

std::vector<unsigned int> itk::simple::ImageRegistrationMethod::m_OptimizerNumberOfSteps
private

Definition at line 748 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerParametersConvergenceTolerance

double itk::simple::ImageRegistrationMethod::m_OptimizerParametersConvergenceTolerance
private

Definition at line 751 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerRelaxationFactor

double itk::simple::ImageRegistrationMethod::m_OptimizerRelaxationFactor
private

Definition at line 737 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerScales

std::vector<double> itk::simple::ImageRegistrationMethod::m_OptimizerScales
private

Definition at line 783 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerScalesCentralRegionRadius

unsigned int itk::simple::ImageRegistrationMethod::m_OptimizerScalesCentralRegionRadius
private

Definition at line 784 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerScalesSmallParameterVariation

double itk::simple::ImageRegistrationMethod::m_OptimizerScalesSmallParameterVariation
private

Definition at line 785 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerScalesType

OptimizerScalesType itk::simple::ImageRegistrationMethod::m_OptimizerScalesType
private

Definition at line 782 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerSeed

unsigned int itk::simple::ImageRegistrationMethod::m_OptimizerSeed
private

Definition at line 761 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerShrinkFactor

double itk::simple::ImageRegistrationMethod::m_OptimizerShrinkFactor
private

Definition at line 760 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerSimplexDelta

double itk::simple::ImageRegistrationMethod::m_OptimizerSimplexDelta
private

Definition at line 750 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerSolutionAccuracy

double itk::simple::ImageRegistrationMethod::m_OptimizerSolutionAccuracy
private

Definition at line 762 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerStepLength

double itk::simple::ImageRegistrationMethod::m_OptimizerStepLength
private

Definition at line 749 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerStepTolerance

double itk::simple::ImageRegistrationMethod::m_OptimizerStepTolerance
private

Definition at line 755 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerTrace

bool itk::simple::ImageRegistrationMethod::m_OptimizerTrace
private

Definition at line 747 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerType

OptimizerType itk::simple::ImageRegistrationMethod::m_OptimizerType
private

Definition at line 727 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerUpperBound

double itk::simple::ImageRegistrationMethod::m_OptimizerUpperBound
private

Definition at line 746 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerValueTolerance

double itk::simple::ImageRegistrationMethod::m_OptimizerValueTolerance
private

Definition at line 756 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerWeights

std::vector<double> itk::simple::ImageRegistrationMethod::m_OptimizerWeights
private

Definition at line 774 of file sitkImageRegistrationMethod.h.

◆ m_OptimizerWithRestarts

bool itk::simple::ImageRegistrationMethod::m_OptimizerWithRestarts
private

Definition at line 753 of file sitkImageRegistrationMethod.h.

◆ m_pfGetCurrentLevel

std::function<unsigned int()> itk::simple::ImageRegistrationMethod::m_pfGetCurrentLevel
private

Definition at line 682 of file sitkImageRegistrationMethod.h.

◆ m_pfGetMetricNumberOfValidPoints

std::function<uint64_t()> itk::simple::ImageRegistrationMethod::m_pfGetMetricNumberOfValidPoints
private

Definition at line 677 of file sitkImageRegistrationMethod.h.

◆ m_pfGetMetricValue

std::function<double()> itk::simple::ImageRegistrationMethod::m_pfGetMetricValue
private

Definition at line 676 of file sitkImageRegistrationMethod.h.

◆ m_pfGetOptimizerConvergenceValue

std::function<double()> itk::simple::ImageRegistrationMethod::m_pfGetOptimizerConvergenceValue
private

Definition at line 675 of file sitkImageRegistrationMethod.h.

◆ m_pfGetOptimizerIteration

std::function<unsigned int()> itk::simple::ImageRegistrationMethod::m_pfGetOptimizerIteration
private

Definition at line 672 of file sitkImageRegistrationMethod.h.

◆ m_pfGetOptimizerLearningRate

std::function<double()> itk::simple::ImageRegistrationMethod::m_pfGetOptimizerLearningRate
private

Definition at line 674 of file sitkImageRegistrationMethod.h.

◆ m_pfGetOptimizerPosition

std::function<std::vector<double>)> itk::simple::ImageRegistrationMethod::m_pfGetOptimizerPosition
private

Definition at line 673 of file sitkImageRegistrationMethod.h.

◆ m_pfGetOptimizerScales

std::function<std::vector<double>)> itk::simple::ImageRegistrationMethod::m_pfGetOptimizerScales
private

Definition at line 678 of file sitkImageRegistrationMethod.h.

◆ m_pfGetOptimizerStopConditionDescription

std::function<std::string()> itk::simple::ImageRegistrationMethod::m_pfGetOptimizerStopConditionDescription
private

Definition at line 679 of file sitkImageRegistrationMethod.h.

◆ m_pfUpdateWithBestValue

std::function<void (itk::TransformBase *outTransform)> itk::simple::ImageRegistrationMethod::m_pfUpdateWithBestValue
private

Definition at line 684 of file sitkImageRegistrationMethod.h.

◆ m_ShrinkFactorsPerLevel

std::vector<unsigned int> itk::simple::ImageRegistrationMethod::m_ShrinkFactorsPerLevel
private

Definition at line 811 of file sitkImageRegistrationMethod.h.

◆ m_SmoothingSigmasAreSpecifiedInPhysicalUnits

bool itk::simple::ImageRegistrationMethod::m_SmoothingSigmasAreSpecifiedInPhysicalUnits
private

Definition at line 813 of file sitkImageRegistrationMethod.h.

◆ m_SmoothingSigmasPerLevel

std::vector<double> itk::simple::ImageRegistrationMethod::m_SmoothingSigmasPerLevel
private

Definition at line 812 of file sitkImageRegistrationMethod.h.

◆ m_StopConditionDescription

std::string itk::simple::ImageRegistrationMethod::m_StopConditionDescription
private

Definition at line 815 of file sitkImageRegistrationMethod.h.

◆ m_TransformBSplineScaleFactors

std::vector<unsigned int> itk::simple::ImageRegistrationMethod::m_TransformBSplineScaleFactors
private

Definition at line 772 of file sitkImageRegistrationMethod.h.

◆ m_VirtualDomainDirection

std::vector<double> itk::simple::ImageRegistrationMethod::m_VirtualDomainDirection
private

Definition at line 713 of file sitkImageRegistrationMethod.h.

◆ m_VirtualDomainOrigin

std::vector<double> itk::simple::ImageRegistrationMethod::m_VirtualDomainOrigin
private

Definition at line 711 of file sitkImageRegistrationMethod.h.

◆ m_VirtualDomainSize

std::vector<uint32_t> itk::simple::ImageRegistrationMethod::m_VirtualDomainSize
private

Definition at line 710 of file sitkImageRegistrationMethod.h.

◆ m_VirtualDomainSpacing

std::vector<double> itk::simple::ImageRegistrationMethod::m_VirtualDomainSpacing
private

Definition at line 712 of file sitkImageRegistrationMethod.h.


The documentation for this class was generated from the following file: