#include <itkComputeDisplacementDistribution.h>
This is a helper class for the automatic parameter estimation of the ASGD optimizer.
More specifically this class computes the Jacobian terms related to the automatic parameter estimation for the adaptive stochastic gradient descent optimizer. Details can be found in the TMI paper
[1] Y. Qiao, B. van Lew, B.P.F. Lelieveldt and M. Staring "Fast Automatic Step Size Estimation for Gradient Descent Optimization of Image Registration," IEEE Transactions on Medical Imaging, vol. 35, no. 2, pp. 391 - 403, February 2016. http://elastix.dev/marius/publications/2016_j_TMIa.php
Definition at line 48 of file itkComputeDisplacementDistribution.h.
Data Structures | |
struct | ComputePerThreadStruct |
struct | MultiThreaderParameterType |
Public Types | |
using | ConstPointer = SmartPointer<const Self> |
using | FixedImageMaskConstPointer = typename FixedImageMaskType::ConstPointer |
using | FixedImageMaskPointer = typename FixedImageMaskType::Pointer |
using | FixedImageMaskType = ImageMaskSpatialObject<Self::FixedImageDimension> |
using | FixedImagePixelType = typename FixedImageType::PixelType |
using | FixedImageRegionType = typename FixedImageType::RegionType |
using | FixedImageType = TFixedImage |
using | NonZeroJacobianIndicesType = typename TransformType::NonZeroJacobianIndicesType |
using | Pointer = SmartPointer<Self> |
using | ScalesType |
using | Self = ComputeDisplacementDistribution |
using | Superclass = ScaledSingleValuedNonLinearOptimizer |
using | TransformPointer = typename TransformType::Pointer |
using | TransformType = TTransform |
Public Types inherited from itk::ScaledSingleValuedNonLinearOptimizer | |
using | ConstPointer = SmartPointer<const Self> |
using | Pointer = SmartPointer<Self> |
using | ScaledCostFunctionPointer = ScaledCostFunctionType::Pointer |
using | ScaledCostFunctionType = ScaledSingleValuedCostFunction |
using | ScalesType = NonLinearOptimizer::ScalesType |
using | Self = ScaledSingleValuedNonLinearOptimizer |
using | Superclass = SingleValuedNonLinearOptimizer |
Public Member Functions | |
virtual void | AfterThreadedCompute (double &jacg, double &maxJJ) |
virtual void | BeforeThreadedCompute (const ParametersType &mu) |
virtual void | Compute (const ParametersType &mu, double &jacg, double &maxJJ, std::string method) |
virtual void | ComputeSingleThreaded (const ParametersType &mu, double &jacg, double &maxJJ, std::string method) |
virtual void | ComputeUsingSearchDirection (const ParametersType &mu, double &jacg, double &maxJJ, std::string methods) |
virtual const char * | GetClassName () const |
virtual const FixedImageMaskType * | GetFixedImageMask () |
virtual const FixedImageRegionType & | GetFixedImageRegion () |
ITK_DISALLOW_COPY_AND_MOVE (ComputeDisplacementDistribution) | |
itkStaticConstMacro (FixedImageDimension, unsigned int, TFixedImage::ImageDimension) | |
virtual void | SetFixedImage (const FixedImageType *_arg) |
virtual void | SetFixedImageMask (const FixedImageMaskType *_arg) |
virtual void | SetFixedImageMask (FixedImageMaskType *_arg) |
void | SetFixedImageRegion (const FixedImageRegionType ®ion) |
virtual void | SetNumberOfJacobianMeasurements (SizeValueType _arg) |
void | SetNumberOfWorkUnits (ThreadIdType numberOfThreads) |
virtual void | SetTransform (TransformType *_arg) |
Public Member Functions inherited from itk::ScaledSingleValuedNonLinearOptimizer | |
const ParametersType & | GetCurrentPosition () const override |
virtual bool | GetMaximize () const |
virtual const ScaledCostFunctionType * | GetScaledCostFunction () |
virtual const ParametersType & | GetScaledCurrentPosition () |
bool | GetUseScales () const |
virtual void | InitializeScales () |
ITK_DISALLOW_COPY_AND_MOVE (ScaledSingleValuedNonLinearOptimizer) | |
virtual void | MaximizeOff () |
virtual void | MaximizeOn () |
void | SetCostFunction (CostFunctionType *costFunction) override |
virtual void | SetMaximize (bool _arg) |
virtual void | SetUseScales (bool arg) |
Static Public Member Functions | |
static Pointer | New () |
Static Public Member Functions inherited from itk::ScaledSingleValuedNonLinearOptimizer | |
static Pointer | New () |
Protected Member Functions | |
ComputeDisplacementDistribution () | |
virtual void | InitializeThreadingParameters () |
itkAlignedTypedef (ITK_CACHE_LINE_ALIGNMENT, PaddedComputePerThreadStruct, AlignedComputePerThreadStruct) | |
itkPadStruct (ITK_CACHE_LINE_ALIGNMENT, ComputePerThreadStruct, PaddedComputePerThreadStruct) | |
void | LaunchComputeThreaderCallback () const |
virtual void | SampleFixedImageForJacobianTerms (ImageSampleContainerPointer &sampleContainer) |
virtual void | ThreadedCompute (ThreadIdType threadID) |
~ComputeDisplacementDistribution () override=default | |
Protected Member Functions inherited from itk::ScaledSingleValuedNonLinearOptimizer | |
virtual void | GetScaledDerivative (const ParametersType ¶meters, DerivativeType &derivative) const |
virtual MeasureType | GetScaledValue (const ParametersType ¶meters) const |
virtual void | GetScaledValueAndDerivative (const ParametersType ¶meters, MeasureType &value, DerivativeType &derivative) const |
void | PrintSelf (std::ostream &os, Indent indent) const override |
ScaledSingleValuedNonLinearOptimizer () | |
void | SetCurrentPosition (const ParametersType ¶m) override |
virtual void | SetScaledCurrentPosition (const ParametersType ¶meters) |
~ScaledSingleValuedNonLinearOptimizer () override=default | |
Static Protected Member Functions | |
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION | ComputeThreaderCallback (void *arg) |
Protected Attributes | |
ScaledSingleValuedCostFunction::Pointer | m_CostFunction {} |
DerivativeType | m_ExactGradient {} |
FixedImageType::ConstPointer | m_FixedImage {} |
FixedImageMaskConstPointer | m_FixedImageMask {} |
FixedImageRegionType | m_FixedImageRegion {} |
SizeValueType | m_NumberOfJacobianMeasurements {} |
SizeValueType | m_NumberOfParameters {} |
MultiThreaderBase::Pointer | m_Threader {} |
TransformPointer | m_Transform {} |
Protected Attributes inherited from itk::ScaledSingleValuedNonLinearOptimizer | |
ScaledCostFunctionPointer | m_ScaledCostFunction {} |
ParametersType | m_ScaledCurrentPosition {} |
Private Attributes | |
std::vector< AlignedComputePerThreadStruct > | m_ComputePerThreadVariables {} |
SizeValueType | m_NumberOfPixelsCounted {} |
ImageSampleContainerPointer | m_SampleContainer {} |
MultiThreaderParameterType | m_ThreaderParameters {} |
bool | m_UseMultiThread {} |
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::ConstPointer = SmartPointer<const Self> |
Definition at line 57 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 176 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 154 of file itkComputeDisplacementDistribution.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::FixedImageMaskConstPointer = typename FixedImageMaskType::ConstPointer |
Definition at line 81 of file itkComputeDisplacementDistribution.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::FixedImageMaskPointer = typename FixedImageMaskType::Pointer |
Definition at line 80 of file itkComputeDisplacementDistribution.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::FixedImageMaskType = ImageMaskSpatialObject<Self::FixedImageDimension> |
Definition at line 79 of file itkComputeDisplacementDistribution.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::FixedImagePixelType = typename FixedImageType::PixelType |
Definition at line 67 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 155 of file itkComputeDisplacementDistribution.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::FixedImageRegionType = typename FixedImageType::RegionType |
Definition at line 70 of file itkComputeDisplacementDistribution.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::FixedImageType = TFixedImage |
typedef
Definition at line 66 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 164 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 163 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 170 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 169 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 167 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 166 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 172 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 171 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 161 of file itkComputeDisplacementDistribution.h.
|
protected |
Samplers.
Definition at line 160 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 156 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 157 of file itkComputeDisplacementDistribution.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::NonZeroJacobianIndicesType = typename TransformType::NonZeroJacobianIndicesType |
Definition at line 82 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 177 of file itkComputeDisplacementDistribution.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::Pointer = SmartPointer<Self> |
Definition at line 56 of file itkComputeDisplacementDistribution.h.
using itk::ScaledSingleValuedNonLinearOptimizer::ScalesType |
Definition at line 85 of file itkScaledSingleValuedNonLinearOptimizer.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::Self = ComputeDisplacementDistribution |
Standard ITK.
Definition at line 54 of file itkComputeDisplacementDistribution.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::Superclass = ScaledSingleValuedNonLinearOptimizer |
Definition at line 55 of file itkComputeDisplacementDistribution.h.
|
protected |
Typedef for multi-threading.
Definition at line 142 of file itkComputeDisplacementDistribution.h.
|
protected |
Typedefs for support of sparse Jacobians and AdvancedTransforms.
Definition at line 175 of file itkComputeDisplacementDistribution.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::TransformPointer = typename TransformType::Pointer |
Definition at line 69 of file itkComputeDisplacementDistribution.h.
using itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::TransformType = TTransform |
Definition at line 68 of file itkComputeDisplacementDistribution.h.
|
protected |
|
overrideprotecteddefault |
|
virtual |
|
virtual |
|
virtual |
The main function that performs the multi-threaded computation.
Reimplemented in itk::ComputePreconditionerUsingDisplacementDistribution< TFixedImage, TTransform >.
|
virtual |
The main function that performs the single-threaded computation.
|
staticprotected |
Compute threader callback function.
|
virtual |
|
virtual |
Run-time type information (and related methods).
Reimplemented from itk::ScaledSingleValuedNonLinearOptimizer.
Reimplemented in itk::ComputePreconditionerUsingDisplacementDistribution< TFixedImage, TTransform >.
|
virtual |
|
virtual |
Get the region over which the metric will be computed.
|
protectedvirtual |
Initialize some multi-threading related parameters.
itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::ITK_DISALLOW_COPY_AND_MOVE | ( | ComputeDisplacementDistribution< TFixedImage, TTransform > | ) |
|
protected |
|
protected |
itk::ComputeDisplacementDistribution< TFixedImage, TTransform >::itkStaticConstMacro | ( | FixedImageDimension | , |
unsigned int | , | ||
TFixedImage::ImageDimension | ) |
Type for the mask of the fixed image. Only pixels that are "inside" this mask will be considered for the computation of the Jacobian terms.
|
protected |
Launch MultiThread Compute.
|
static |
Method for creation through the object factory.
|
protectedvirtual |
Sample the fixed image to compute the Jacobian terms.
|
virtual |
Set the fixed image.
|
virtual |
|
virtual |
Set/Get the fixed image mask.
|
inline |
Set the region over which the metric will be computed.
Definition at line 100 of file itkComputeDisplacementDistribution.h.
|
virtual |
Set some parameters.
|
inline |
Set the number of threads.
Definition at line 125 of file itkComputeDisplacementDistribution.h.
|
virtual |
Set the transform.
|
protectedvirtual |
The threaded implementation of Compute().
|
mutableprivate |
Definition at line 220 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 148 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 150 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 144 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 146 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 145 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 149 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 151 of file itkComputeDisplacementDistribution.h.
|
private |
Definition at line 222 of file itkComputeDisplacementDistribution.h.
|
private |
Definition at line 224 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 152 of file itkComputeDisplacementDistribution.h.
|
mutableprivate |
Definition at line 218 of file itkComputeDisplacementDistribution.h.
|
protected |
Definition at line 147 of file itkComputeDisplacementDistribution.h.
|
private |
Definition at line 223 of file itkComputeDisplacementDistribution.h.
Generated on 2024-07-17 for elastix by 1.11.0 (9b424b03c9833626cd435af22a444888fbbb192d) |