go home Home | Main Page | Topics | Namespace List | Class Hierarchy | Alphabetical List | Data Structures | File List | Namespace Members | Data Fields | Globals | Related Pages
itk::ImpactImageToImageMetric< TFixedImage, TMovingImage > Class Template Reference

#include <itkImpactImageToImageMetric.h>

Detailed Description

template<typename TFixedImage, typename TMovingImage>
class itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >

A semantic similarity metric for multimodal image registration based on deep learning features.

This class define a loss by compares the fixed and moving images using high-level semantic representations extracted from pretrained deep learning models, rather than relying on raw pixel intensities or handcrafted features.

Author
V. Boussot, Univ. Rennes, INSERM, LTSI- UMR 1099, F-35000 Rennes, France
Note
This work was funded by the French National Research Agency as part of the VATSop project (ANR-20-CE19-0015).
If you use the Impact anywhere we would appreciate if you cite the following article:
V. Boussot et al., IMPACT: A Generic Semantic Loss for Multimodal Medical Image Registration, arXiv preprint arXiv:2503.24121 (2025). https://doi.org/10.48550/arXiv.2503.24121

Definition at line 58 of file itkImpactImageToImageMetric.h.

Inheritance diagram for itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >:

Classes

struct  FeaturesMaps
struct  LossPerThreadStruct

Public Types

using ConstPointer = SmartPointer<const Self>
using DerivativeValueType
using FixedImageLimiterOutputType
using FixedImageLimiterType
using FixedImageMaskPointer
using FixedImageMaskType
using FixedImagePixelType
using ImageSampleContainerPointer
using ImageSampleContainerType
using ImageSamplerPointer
using ImageSamplerType
using MovingImageDerivativeScalesType
using MovingImageLimiterOutputType
using MovingImageLimiterType
using MovingImageMaskPointer
using MovingImageMaskType
using MovingImageRegionType
using NumberOfParametersType
using Pointer = SmartPointer<Self>
using Self = ImpactImageToImageMetric
using Superclass = AdvancedImageToImageMetric<TFixedImage, TMovingImage>
using ThreadInfoType
Public Types inherited from itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >
using AdvancedTransformType = AdvancedTransform<ScalarType, FixedImageDimension, MovingImageDimension>
using BSplineOrder1TransformPointer = typename BSplineOrder1TransformType::Pointer
using BSplineOrder1TransformType = AdvancedBSplineDeformableTransform<ScalarType, FixedImageDimension, 1>
using BSplineOrder2TransformPointer = typename BSplineOrder2TransformType::Pointer
using BSplineOrder2TransformType = AdvancedBSplineDeformableTransform<ScalarType, FixedImageDimension, 2>
using BSplineOrder3TransformPointer = typename BSplineOrder3TransformType::Pointer
using BSplineOrder3TransformType = AdvancedBSplineDeformableTransform<ScalarType, FixedImageDimension, 3>
using CombinationTransformType = AdvancedCombinationTransform<ScalarType, FixedImageDimension>
using ConstPointer = SmartPointer<const Self>
using DerivativeValueType = typename DerivativeType::ValueType
using FixedImageLimiterOutputType = typename FixedImageLimiterType::OutputType
using FixedImageLimiterPointer = typename FixedImageLimiterType::Pointer
using FixedImageLimiterType = LimiterFunctionBase<RealType, FixedImageDimension>
using FixedImageMaskConstPointer = SmartPointer<const FixedImageMaskType>
using FixedImageMaskPointer = SmartPointer<FixedImageMaskType>
using FixedImageMaskType = ImageMaskSpatialObject<Self::FixedImageDimension>
using FixedImagePixelType = typename FixedImageType::PixelType
using FixedImagePointer = typename FixedImageType::Pointer
using ImageSampleContainerPointer = typename ImageSamplerType::OutputVectorContainerPointer
using ImageSampleContainerType = typename ImageSamplerType::OutputVectorContainerType
using ImageSamplerPointer = typename ImageSamplerType::Pointer
using ImageSamplerType = ImageSamplerBase<FixedImageType>
using MovingImageDerivativeScalesType = FixedArray<double, Self::MovingImageDimension>
using MovingImageLimiterOutputType = typename MovingImageLimiterType::OutputType
using MovingImageLimiterPointer = typename MovingImageLimiterType::Pointer
using MovingImageLimiterType = LimiterFunctionBase<RealType, MovingImageDimension>
using MovingImageMaskConstPointer = SmartPointer<const MovingImageMaskType>
using MovingImageMaskPointer = SmartPointer<MovingImageMaskType>
using MovingImageMaskType = ImageMaskSpatialObject<Self::MovingImageDimension>
using MovingImagePointer = typename MovingImageType::Pointer
using MovingImageRegionType = typename MovingImageType::RegionType
using NumberOfParametersType = typename AdvancedTransformType::NumberOfParametersType
using Pointer = SmartPointer<Self>
using ScalarType = typename TransformType::ScalarType
using Self = AdvancedImageToImageMetric
using Superclass = ImageToImageMetric<TFixedImage, TMovingImage>
using ThreadInfoType = MultiThreaderBase::WorkUnitInfo

Public Member Functions

virtual const char * GetClassName () const
virtual unsigned int GetCurrentLevel () const
void GetDerivative (const ParametersType &parameters, DerivativeType &derivative) const override
virtual torch::Device GetDevice () const
virtual std::vector< std::string > GetDistance () const
virtual std::string GetFeatureMapsPath () const
virtual int GetFeaturesMapUpdateInterval () const
virtual const std::vector< ImpactModelConfiguration > & GetFixedModelsConfiguration ()
virtual std::vector< floatGetLayersWeight () const
virtual std::string GetMode () const
virtual const std::vector< ImpactModelConfiguration > & GetMovingModelsConfiguration ()
virtual std::vector< unsigned intGetPCA () const
virtual unsigned int GetSeed () const
virtual std::vector< unsigned intGetSubsetFeatures () const
virtual bool GetUseMixedPrecision () const
MeasureType GetValue (const ParametersType &parameters) const override
void GetValueAndDerivative (const ParametersType &parameters, MeasureType &value, DerivativeType &derivative) const override
void GetValueAndDerivativeSingleThreaded (const ParametersType &parameters, MeasureType &value, DerivativeType &derivative) const
virtual MeasureType GetValueSingleThreaded (const ParametersType &parameters) const
virtual bool GetWriteFeatureMaps () const
void Initialize () override
 ITK_DISALLOW_COPY_AND_MOVE (ImpactImageToImageMetric)
 itkStaticConstMacro (FixedImageDimension, unsigned int, FixedImageType::ImageDimension)
 itkStaticConstMacro (MovingImageDimension, unsigned int, MovingImageType::ImageDimension)
virtual void SetCurrentLevel (unsigned int _arg)
virtual void SetDevice (torch::Device _arg)
virtual void SetDistance (std::vector< std::string > _arg)
virtual void SetFeatureMapsPath (std::string _arg)
virtual void SetFeaturesMapUpdateInterval (int _arg)
virtual void SetFixedModelsConfiguration (std::vector< ImpactModelConfiguration > _arg)
virtual void SetLayersWeight (std::vector< float > _arg)
virtual void SetMode (std::string _arg)
virtual void SetMovingModelsConfiguration (std::vector< ImpactModelConfiguration > _arg)
virtual void SetPCA (std::vector< unsigned int > _arg)
virtual void SetSeed (unsigned int _arg)
virtual void SetSubsetFeatures (std::vector< unsigned int > _arg)
virtual void SetUseMixedPrecision (bool _arg)
virtual void SetWriteFeatureMaps (bool _arg)
Public Member Functions inherited from itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >
virtual void BeforeThreadedGetValueAndDerivative (const TransformParametersType &parameters) const
virtual const FixedImageLimiterTypeGetFixedImageLimiter ()
const FixedImageMaskTypeGetFixedImageMask () const override
virtual double GetFixedLimitRangeRatio () const
ImageSamplerTypeGetImageSampler () const
virtual const MovingImageDerivativeScalesTypeGetMovingImageDerivativeScales ()
virtual const MovingImageLimiterTypeGetMovingImageLimiter ()
const MovingImageMaskTypeGetMovingImageMask () const override
virtual double GetMovingLimitRangeRatio () const
virtual double GetRequiredRatioOfValidSamples () const
virtual bool GetScaleGradientWithRespectToMovingImageOrientation () const
const AdvancedTransformTypeGetTransform () const override
AdvancedTransformTypeGetTransform () override
virtual bool GetUseFixedImageLimiter () const
virtual bool GetUseImageSampler () const
virtual const boolGetUseMetricSingleThreaded ()
virtual bool GetUseMovingImageDerivativeScales () const
virtual bool GetUseMovingImageLimiter () const
virtual const boolGetUseMultiThread ()
void Initialize () override
 ITK_DISALLOW_COPY_AND_MOVE (AdvancedImageToImageMetric)
 itkOverrideGetNameOfClassMacro (AdvancedImageToImageMetric)
 itkStaticConstMacro (FixedImageDimension, unsigned int, TFixedImage::ImageDimension)
 itkStaticConstMacro (MovingImageDimension, unsigned int, TMovingImage::ImageDimension)
virtual void SetFixedImageLimiter (FixedImageLimiterType *_arg)
virtual void SetFixedImageMask (const FixedImageMaskType *const arg)
virtual void SetFixedLimitRangeRatio (double _arg)
virtual void SetImageSampler (ImageSamplerType *_arg)
virtual void SetMovingImageDerivativeScales (MovingImageDerivativeScalesType _arg)
virtual void SetMovingImageLimiter (MovingImageLimiterType *_arg)
virtual void SetMovingImageMask (const MovingImageMaskType *const arg)
virtual void SetMovingLimitRangeRatio (double _arg)
void SetRandomVariateGenerator (Statistics::MersenneTwisterRandomVariateGenerator &randomVariateGenerator)
virtual void SetRequiredRatioOfValidSamples (double _arg)
virtual void SetScaleGradientWithRespectToMovingImageOrientation (bool _arg)
virtual void SetTransform (AdvancedTransformType *arg)
virtual void SetUseMetricSingleThreaded (bool _arg)
virtual void SetUseMovingImageDerivativeScales (bool _arg)
virtual void SetUseMultiThread (bool _arg)
virtual void UseMetricSingleThreadedOff ()
virtual void UseMetricSingleThreadedOn ()
virtual void UseMultiThreadOff ()
virtual void UseMultiThreadOn ()

Static Public Member Functions

static Pointer New ()

Protected Types

using BSplineInterpolatorType
using FixedImageIndexType
using FixedImageIndexValueType
using FixedImagePointType
using MovingImageContinuousIndexType
using MovingImageDerivativeType
using MovingImageIndexType
using MovingImagePointType
using NonZeroJacobianIndicesType
Protected Types inherited from itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >
using BSplineInterpolatorFloatPointer = typename BSplineInterpolatorFloatType::Pointer
using BSplineInterpolatorFloatType
using BSplineInterpolatorPointer = typename BSplineInterpolatorType::Pointer
using BSplineInterpolatorType
using FixedImageIndexType = typename FixedImageType::IndexType
using FixedImageIndexValueType = typename FixedImageIndexType::IndexValueType
using FixedImagePointType = typename TransformType::InputPointType
using LinearInterpolatorPointer = typename LinearInterpolatorType::Pointer
using LinearInterpolatorType = AdvancedLinearInterpolateImageFunction<MovingImageType, CoordinateRepresentationType>
using MovingImageContinuousIndexType = typename InterpolatorType::ContinuousIndexType
using MovingImageDerivativeType = typename BSplineInterpolatorType::CovariantVectorType
using MovingImageIndexType = typename MovingImageType::IndexType
using MovingImagePointType = typename TransformType::OutputPointType
using NonZeroJacobianIndicesType = typename AdvancedTransformType::NonZeroJacobianIndicesType
using ReducedBSplineInterpolatorPointer = typename ReducedBSplineInterpolatorType::Pointer
using ReducedBSplineInterpolatorType

Protected Member Functions

void AfterThreadedGetValue (MeasureType &value) const override
void AfterThreadedGetValueAndDerivative (MeasureType &value, DerivativeType &derivative) const override
unsigned int ComputeValue (const std::vector< FixedImagePointType > &fixedPoints, LossPerThreadStruct &loss) const
unsigned int ComputeValueAndDerivativeJacobian (const std::vector< FixedImagePointType > &fixedPoints, LossPerThreadStruct &loss) const
unsigned int ComputeValueAndDerivativeStatic (const std::vector< FixedImagePointType > &fixedPoints, LossPerThreadStruct &loss) const
unsigned int ComputeValueStatic (const std::vector< FixedImagePointType > &fixedPoints, LossPerThreadStruct &loss) const
 ImpactImageToImageMetric ()
void InitializeThreadingParameters () const override
bool SampleCheck (const FixedImagePointType &fixedImageCenterCoordinate) const
bool SampleCheck (const FixedImagePointType &fixedImageCenterCoordinate, const std::vector< std::vector< float > > &patchIndex) const
void ThreadedGetValue (ThreadIdType threadID) const override
void ThreadedGetValueAndDerivative (ThreadIdType threadID) const override
void UpdateFixedFeaturesMaps ()
void UpdateMovingFeaturesMaps ()
 ~ImpactImageToImageMetric () override=default
Protected Member Functions inherited from itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >
 AdvancedImageToImageMetric ()
void CheckForAdvancedTransform ()
void CheckForBSplineInterpolator ()
void CheckForBSplineTransform () const
void CheckNumberOfSamples () const
virtual bool EvaluateMovingImageValueAndDerivative (const MovingImagePointType &mappedPoint, RealType &movingImageValue, MovingImageDerivativeType *gradient) const
bool EvaluateTransformJacobian (const FixedImagePointType &fixedImagePoint, TransformJacobianType &jacobian, NonZeroJacobianIndicesType &nzji) const
virtual void EvaluateTransformJacobianInnerProduct (const TransformJacobianType &jacobian, const MovingImageDerivativeType &movingImageDerivative, DerivativeType &imageJacobian) const
bool FastEvaluateMovingImageValueAndDerivative (const MovingImagePointType &mappedPoint, RealType &movingImageValue, MovingImageDerivativeType *gradient, const ThreadIdType threadId) const
Statistics::MersenneTwisterRandomVariateGenerator & GetMutableRandomVariateGenerator () const
Statistics::MersenneTwisterRandomVariateGenerator & GetRandomVariateGenerator ()
virtual void InitializeImageSampler ()
void InitializeLimiters ()
virtual bool IsInsideMovingMask (const MovingImagePointType &point) const
 itkAlignedTypedef (ITK_CACHE_LINE_ALIGNMENT, PaddedGetValueAndDerivativePerThreadStruct, AlignedGetValueAndDerivativePerThreadStruct)
 itkPadStruct (ITK_CACHE_LINE_ALIGNMENT, GetValueAndDerivativePerThreadStruct, PaddedGetValueAndDerivativePerThreadStruct)
void LaunchGetValueAndDerivativeThreaderCallback () const
void LaunchGetValueThreaderCallback () const
void PrintSelf (std::ostream &os, Indent indent) const override
void SetFixedImageMask (const typename Superclass::FixedImageMaskType *) final
void SetFixedImageMask (typename Superclass::FixedImageMaskType *) final
void SetMovingImageMask (const typename Superclass::MovingImageMaskType *) final
void SetMovingImageMask (typename Superclass::MovingImageMaskType *) final
virtual void SetUseFixedImageLimiter (bool _arg)
virtual void SetUseImageSampler (bool _arg)
virtual void SetUseMovingImageLimiter (bool _arg)
MovingImagePointType TransformPoint (const FixedImagePointType &fixedImagePoint) const
 ~AdvancedImageToImageMetric () override=default

Private Types

using FeaturesImageType = itk::VectorImage<float, FixedImageDimension>
using FeaturesInterpolatorType
using FeaturesMaps = typename ImpactImageToImageMetric<TFixedImage, TMovingImage>::FeaturesMaps
using FixedInterpolatorType = BSplineInterpolateImageFunction<FixedImageType, CoordinateRepresentationType, double>

Private Member Functions

torch::Tensor EvaluateFixedImagesPatchValue (const FixedImagePointType &fixedImageCenterCoordinate, const std::vector< std::vector< float > > &patchIndex, const std::vector< int64_t > &patchSize) const
torch::Tensor EvaluateMovingImagesPatchValue (const FixedImagePointType &fixedImageCenterCoordinate, const std::vector< std::vector< float > > &patchIndex, const std::vector< int64_t > &patchSize) const
torch::Tensor EvaluateMovingImagesPatchValuesAndJacobians (const FixedImagePointType &fixedImageCenterCoordinate, torch::Tensor &movingImagesPatchesJacobians, const std::vector< std::vector< float > > &patchIndex, const std::vector< int64_t > &patchSize, int s) const
template<typename ImagePointType>
std::vector< ImagePointType > GeneratePatchIndex (const std::vector< ImpactModelConfiguration > &modelConfig, std::mt19937 &randomGenerator, const std::vector< ImagePointType > &fixedPointsTmp, std::vector< std::vector< std::vector< std::vector< float > > > > &patchIndex) const
std::vector< unsigned intGetSubsetOfFeatures (const std::vector< unsigned int > &featuresIndex, std::mt19937 &randomGenerator, int n) const
 itkAlignedTypedef (ITK_CACHE_LINE_ALIGNMENT, PaddedLossPerThreadStruct, AlignedLossPerThreadStruct)
 itkPadStruct (ITK_CACHE_LINE_ALIGNMENT, LossPerThreadStruct, PaddedLossPerThreadStruct)

Private Attributes

unsigned int m_CurrentLevel
torch::Device m_Device = torch::Device(torch::kCPU)
std::vector< std::string > m_Distance
std::string m_FeatureMapsPath
std::vector< std::vector< unsigned int > > m_FeaturesIndexes
int m_FeaturesMapUpdateInterval
std::vector< FeaturesMapsm_FixedFeaturesMaps
InterpolatorPointer m_FixedInterpolator
std::vector< ImpactModelConfigurationm_FixedModelsConfiguration
std::vector< floatm_LayersWeight
std::unique_ptr< AlignedLossPerThreadStruct[]> m_LossThreadStruct { nullptr }
int m_LossThreadStructSize = 0
std::string m_Mode
std::vector< FeaturesMapsm_MovingFeaturesMaps
std::vector< ImpactModelConfigurationm_MovingModelsConfiguration
std::vector< unsigned intm_PCA
std::vector< torch::Tensor > m_PrincipalComponents
unsigned int m_Seed
std::vector< unsigned intm_SubsetFeatures
bool m_UseMixedPrecision
bool m_WriteFeatureMaps

Additional Inherited Members

Static Protected Member Functions inherited from itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION AccumulateDerivativesThreaderCallback (void *arg)
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION GetValueAndDerivativeThreaderCallback (void *arg)
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION GetValueThreaderCallback (void *arg)
Protected Attributes inherited from itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >
AdvancedTransformType::Pointer m_AdvancedTransform { nullptr }
FixedImageLimiterOutputType m_FixedImageMaxLimit { 1 }
FixedImageLimiterOutputType m_FixedImageMinLimit { 0 }
FixedImagePixelType m_FixedImageTrueMax { 1 }
FixedImagePixelType m_FixedImageTrueMin { 0 }
double m_FixedLimitRangeRatio { 0.01 }
std::unique_ptr< AlignedGetValueAndDerivativePerThreadStruct[]> m_GetValueAndDerivativePerThreadVariables
ThreadIdType m_GetValueAndDerivativePerThreadVariablesSize { 0 }
ImageSamplerPointer m_ImageSampler { nullptr }
MovingImageLimiterOutputType m_MovingImageMaxLimit { 1 }
MovingImageLimiterOutputType m_MovingImageMinLimit { 0 }
MovingImagePixelType m_MovingImageTrueMax { 1 }
MovingImagePixelType m_MovingImageTrueMin { 0 }
double m_MovingLimitRangeRatio { 0.01 }
MultiThreaderParameterType m_ThreaderMetricParameters {}
bool m_TransformIsBSpline { false }
bool m_UseMetricSingleThreaded { true }
bool m_UseMultiThread { false }

Member Typedef Documentation

◆ BSplineInterpolatorType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::BSplineInterpolatorType
protected

Typedefs used for computing image derivatives.

Definition at line 350 of file itkAdvancedImageToImageMetric.h.

◆ ConstPointer

template<typename TFixedImage, typename TMovingImage>
using itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::ConstPointer = SmartPointer<const Self>

Definition at line 67 of file itkImpactImageToImageMetric.h.

◆ DerivativeValueType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::DerivativeValueType

Definition at line 135 of file itkAdvancedImageToImageMetric.h.

◆ FeaturesImageType

template<typename TFixedImage, typename TMovingImage>
using itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::FeaturesImageType = itk::VectorImage<float, FixedImageDimension>
private

Feature maps are stored as VectorImages of floats with same dimension as fixed image.

Definition at line 537 of file itkImpactImageToImageMetric.h.

◆ FeaturesInterpolatorType

template<typename TFixedImage, typename TMovingImage>
using itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::FeaturesInterpolatorType
private
Initial value:
BSplineInterpolateImageFunction<itk::Image<float, FixedImageDimension>, CoordinateRepresentationType, float>>
Helper class to interpolate each component of a VectorImage using separate B-Spline interpolators.
itk::VectorImage< float, FixedImageDimension > FeaturesImageType

Interpolator for feature maps (vector-valued), using scalar B-spline interpolation.

Definition at line 539 of file itkImpactImageToImageMetric.h.

◆ FeaturesMaps

template<typename TFixedImage, typename TMovingImage>
using itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::FeaturesMaps = typename ImpactImageToImageMetric<TFixedImage, TMovingImage>::FeaturesMaps
private

Definition at line 564 of file itkImpactImageToImageMetric.h.

◆ FixedImageIndexType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::FixedImageIndexType
protected

Protected Typedefs Typedefs for indices and points.

Definition at line 342 of file itkAdvancedImageToImageMetric.h.

◆ FixedImageIndexValueType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::FixedImageIndexValueType
protected

Definition at line 343 of file itkAdvancedImageToImageMetric.h.

◆ FixedImageLimiterOutputType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::FixedImageLimiterOutputType

Definition at line 152 of file itkAdvancedImageToImageMetric.h.

◆ FixedImageLimiterType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::FixedImageLimiterType

Typedefs for Limiter support.

Definition at line 150 of file itkAdvancedImageToImageMetric.h.

◆ FixedImageMaskPointer

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::FixedImageMaskPointer

Definition at line 127 of file itkAdvancedImageToImageMetric.h.

◆ FixedImageMaskType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::FixedImageMaskType

Definition at line 126 of file itkAdvancedImageToImageMetric.h.

◆ FixedImagePixelType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::FixedImagePixelType

Some useful extra typedefs.

Definition at line 139 of file itkAdvancedImageToImageMetric.h.

◆ FixedImagePointType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::FixedImagePointType
protected

Definition at line 345 of file itkAdvancedImageToImageMetric.h.

◆ FixedInterpolatorType

template<typename TFixedImage, typename TMovingImage>
using itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::FixedInterpolatorType = BSplineInterpolateImageFunction<FixedImageType, CoordinateRepresentationType, double>
private

Interpolator for fixed image intensities, using B-spline of order 3 (double precision).

Definition at line 535 of file itkImpactImageToImageMetric.h.

◆ ImageSampleContainerPointer

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::ImageSampleContainerPointer

Definition at line 147 of file itkAdvancedImageToImageMetric.h.

◆ ImageSampleContainerType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::ImageSampleContainerType

Definition at line 146 of file itkAdvancedImageToImageMetric.h.

◆ ImageSamplerPointer

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::ImageSamplerPointer

Definition at line 145 of file itkAdvancedImageToImageMetric.h.

◆ ImageSamplerType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::ImageSamplerType

Typedefs for the ImageSampler.

Definition at line 144 of file itkAdvancedImageToImageMetric.h.

◆ MovingImageContinuousIndexType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::MovingImageContinuousIndexType
protected

Definition at line 347 of file itkAdvancedImageToImageMetric.h.

◆ MovingImageDerivativeScalesType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::MovingImageDerivativeScalesType

Definition at line 141 of file itkAdvancedImageToImageMetric.h.

◆ MovingImageDerivativeType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::MovingImageDerivativeType
protected

Definition at line 361 of file itkAdvancedImageToImageMetric.h.

◆ MovingImageIndexType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::MovingImageIndexType
protected

Definition at line 344 of file itkAdvancedImageToImageMetric.h.

◆ MovingImageLimiterOutputType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::MovingImageLimiterOutputType

Definition at line 155 of file itkAdvancedImageToImageMetric.h.

◆ MovingImageLimiterType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::MovingImageLimiterType

Definition at line 153 of file itkAdvancedImageToImageMetric.h.

◆ MovingImageMaskPointer

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::MovingImageMaskPointer

Definition at line 130 of file itkAdvancedImageToImageMetric.h.

◆ MovingImageMaskType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::MovingImageMaskType

Definition at line 129 of file itkAdvancedImageToImageMetric.h.

◆ MovingImagePointType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::MovingImagePointType
protected

Definition at line 346 of file itkAdvancedImageToImageMetric.h.

◆ MovingImageRegionType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::MovingImageRegionType

Definition at line 140 of file itkAdvancedImageToImageMetric.h.

◆ NonZeroJacobianIndicesType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::NonZeroJacobianIndicesType
protected

Typedefs for support of sparse Jacobians and compact support of transformations.

Definition at line 364 of file itkAdvancedImageToImageMetric.h.

◆ NumberOfParametersType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::NumberOfParametersType

Definition at line 160 of file itkAdvancedImageToImageMetric.h.

◆ Pointer

template<typename TFixedImage, typename TMovingImage>
using itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::Pointer = SmartPointer<Self>

Definition at line 66 of file itkImpactImageToImageMetric.h.

◆ Self

template<typename TFixedImage, typename TMovingImage>
using itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::Self = ImpactImageToImageMetric

Standard class typedefs.

Definition at line 64 of file itkImpactImageToImageMetric.h.

◆ Superclass

template<typename TFixedImage, typename TMovingImage>
using itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::Superclass = AdvancedImageToImageMetric<TFixedImage, TMovingImage>

Definition at line 65 of file itkImpactImageToImageMetric.h.

◆ ThreadInfoType

template<typename TFixedImage, typename TMovingImage>
using itk::AdvancedImageToImageMetric< TFixedImage, TMovingImage >::ThreadInfoType

Typedef for multi-threading.

Definition at line 172 of file itkAdvancedImageToImageMetric.h.

Constructor & Destructor Documentation

◆ ImpactImageToImageMetric()

template<typename TFixedImage, typename TMovingImage>
itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::ImpactImageToImageMetric ( )
protected

◆ ~ImpactImageToImageMetric()

template<typename TFixedImage, typename TMovingImage>
itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::~ImpactImageToImageMetric ( )
overrideprotecteddefault

Member Function Documentation

◆ AfterThreadedGetValue()

template<typename TFixedImage, typename TMovingImage>
void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::AfterThreadedGetValue ( MeasureType & value) const
overrideprotectedvirtual

Combines the similarity values computed by all threads.

This method aggregates the loss contributions stored in each thread’s LossPerThreadStruct and combines them into a global scalar value. This aggregated value is then used by the optimizer during the optimization process.

Parameters
valueThe scalar value to store the aggregated similarity result.

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

◆ AfterThreadedGetValueAndDerivative()

template<typename TFixedImage, typename TMovingImage>
void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::AfterThreadedGetValueAndDerivative ( MeasureType & value,
DerivativeType & derivative ) const
overrideprotectedvirtual

Combines the values and gradients computed by all threads.

This method performs the final reduction step to aggregate the values and gradients computed by each thread into global loss and gradient vectors. These aggregated results are then used in the optimization process.

Parameters
valueThe global loss value computed by combining the contributions from all threads.
derivativeThe global gradient vector computed by combining the gradients from all threads.

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

◆ ComputeValue()

template<typename TFixedImage, typename TMovingImage>
unsigned int itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::ComputeValue ( const std::vector< FixedImagePointType > & fixedPoints,
LossPerThreadStruct & loss ) const
protected

Computes the semantic similarity value using the current transform parameters.

This method computes the loss at all sampled points using the current transformation parameters, without calculating derivatives. It performs patch-based inference and feature comparison in Jacobian mode.

Parameters
fixedPointsA vector of sampled points in the fixed image.
lossThe loss objects (one per layer) where the values will be accumulated.
Returns
The number of valid samples used in the similarity computation.

◆ ComputeValueAndDerivativeJacobian()

template<typename TFixedImage, typename TMovingImage>
unsigned int itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::ComputeValueAndDerivativeJacobian ( const std::vector< FixedImagePointType > & fixedPoints,
LossPerThreadStruct & loss ) const
protected

Computes both the semantic similarity value and its derivative using Jacobian mode.

This method computes both the similarity value and its derivative (gradient) by backpropagating through the model. This allows the metric to assess the sensitivity to transformation parameters, which is crucial for gradient-based optimization.

Parameters
fixedPointsA vector of sampled points in the fixed image.
lossLoss objects to store both the values and gradients for each semantic layer.
Returns
The number of valid samples used in the similarity and gradient computation.

◆ ComputeValueAndDerivativeStatic()

template<typename TFixedImage, typename TMovingImage>
unsigned int itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::ComputeValueAndDerivativeStatic ( const std::vector< FixedImagePointType > & fixedPoints,
LossPerThreadStruct & loss ) const
protected

Computes the value and derivative in static mode (using precomputed features).

This method computes the similarity value and its derivative (gradient) using precomputed feature maps. Gradients are computed via the chain rule, using the interpolated feature fields rather than backpropagating through the model.

Parameters
fixedPointsA vector of sampled points in the fixed image.
lossLoss objects to store both the values and gradients for each semantic layer.
Returns
The number of valid samples used in the similarity and gradient computation.

◆ ComputeValueStatic()

template<typename TFixedImage, typename TMovingImage>
unsigned int itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::ComputeValueStatic ( const std::vector< FixedImagePointType > & fixedPoints,
LossPerThreadStruct & loss ) const
protected

Computes the semantic similarity value in static mode (using precomputed feature maps).

Unlike ComputeValue(), this method uses pre-extracted static feature maps for both the fixed and moving images, avoiding repeated forward passes through the model

Parameters
fixedPointsA vector of sampled points in the fixed image.
lossThe loss objects (one per semantic layer) where the values will be accumulated.
Returns
The number of valid samples used in the similarity computation.

◆ EvaluateFixedImagesPatchValue()

template<typename TFixedImage, typename TMovingImage>
torch::Tensor itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::EvaluateFixedImagesPatchValue ( const FixedImagePointType & fixedImageCenterCoordinate,
const std::vector< std::vector< float > > & patchIndex,
const std::vector< int64_t > & patchSize ) const
private

Extracts a fixed image patch tensor centered at a given point using the precomputed patch index.

This method extracts a patch from the fixed image centered at the specified point, based on the provided patchIndex and patchSize. Interpolation is performed using the fixed image interpolator to sample the feature values at the specified coordinates.

Parameters
fixedImageCenterCoordinateThe coordinate of the center point of the patch in the fixed image.
patchIndexA vector defining the relative indices for extracting the patch.
patchSizeThe size of the patch to extract.
Returns
A tensor representing the extracted patch from the fixed image.

◆ EvaluateMovingImagesPatchValue()

template<typename TFixedImage, typename TMovingImage>
torch::Tensor itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::EvaluateMovingImagesPatchValue ( const FixedImagePointType & fixedImageCenterCoordinate,
const std::vector< std::vector< float > > & patchIndex,
const std::vector< int64_t > & patchSize ) const
private

Extracts a moving image patch tensor (intensity values) corresponding to a fixed point.

This method extracts a patch from the moving image, centered at the corresponding transformed point from the fixed image. The extraction is performed using the provided transform and the moving image interpolator to sample intensity values at the transformed coordinates.

Parameters
fixedImageCenterCoordinateThe coordinate of the center point in the fixed image.
patchIndexA 2D vector defining the relative indices for extracting the patch.
patchSizeThe size of the patch to extract.
Returns
A tensor representing the extracted patch from the moving image.

◆ EvaluateMovingImagesPatchValuesAndJacobians()

template<typename TFixedImage, typename TMovingImage>
torch::Tensor itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::EvaluateMovingImagesPatchValuesAndJacobians ( const FixedImagePointType & fixedImageCenterCoordinate,
torch::Tensor & movingImagesPatchesJacobians,
const std::vector< std::vector< float > > & patchIndex,
const std::vector< int64_t > & patchSize,
int s ) const
private

Extracts moving image patch values and computes the spatial Jacobians with respect to image coordinates.

This method extracts a patch from the moving image centered at the corresponding transformed point from the fixed image. It also computes the image gradient. This is used in Jacobian mode to enable gradient computation via the chain rule during optimization.

Parameters
fixedImageCenterCoordinateThe coordinate of the center point in the fixed image.
movingImagesPatchesJacobiansA tensor to store the computed Jacobians (image gradients) for each patch.
patchIndexA vector defining the relative indices for extracting the patch.
patchSizeThe size of the patch to extract.
sAn index to determine where to store the computed Jacobians in movingImagesPatchesJacobians.
Returns
A tensor representing the extracted patch from the moving image.

◆ GeneratePatchIndex()

template<typename TFixedImage, typename TMovingImage>
template<typename ImagePointType>
std::vector< ImagePointType > itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GeneratePatchIndex ( const std::vector< ImpactModelConfiguration > & modelConfig,
std::mt19937 & randomGenerator,
const std::vector< ImagePointType > & fixedPointsTmp,
std::vector< std::vector< std::vector< std::vector< float > > > > & patchIndex ) const
private

Generates valid patch indices and filters out invalid points.

This method takes a list of fixed points and model configurations to generate patch indices for each point. It filters out points that lie outside the mask or image boundaries, returning only the valid fixed points. The valid patch indices are stored in patchIndex.

Parameters
modelConfigA vector of model configurations, each specifying the properties of the models used.
randomGeneratorA random number generator used for random sampling in the patch generation process.
fixedPointsTmpA vector of fixed points to generate patch indices for.
patchIndexA reference to a 4D vector to store the generated patch indices for each valid fixed point.
Returns
A vector of valid fixed points that are inside the mask or boundaries of the image.

◆ GetClassName()

template<typename TFixedImage, typename TMovingImage>
virtual const char * itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetClassName ( ) const
virtual

Run-time type information (and related methods).

◆ GetCurrentLevel()

template<typename TFixedImage, typename TMovingImage>
virtual unsigned int itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetCurrentLevel ( ) const
virtual

◆ GetDerivative()

template<typename TFixedImage, typename TMovingImage>
void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetDerivative ( const ParametersType & parameters,
DerivativeType & derivative ) const
override

Compute the gradient (derivative) of the similarity value with respect to transformation parameters. Used in gradient-based optimization methods. Internally supports multi-threaded computation.

◆ GetDevice()

template<typename TFixedImage, typename TMovingImage>
virtual torch::Device itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetDevice ( ) const
virtual

◆ GetDistance()

template<typename TFixedImage, typename TMovingImage>
virtual std::vector< std::string > itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetDistance ( ) const
virtual

◆ GetFeatureMapsPath()

template<typename TFixedImage, typename TMovingImage>
virtual std::string itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetFeatureMapsPath ( ) const
virtual

◆ GetFeaturesMapUpdateInterval()

template<typename TFixedImage, typename TMovingImage>
virtual int itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetFeaturesMapUpdateInterval ( ) const
virtual

◆ GetFixedModelsConfiguration()

template<typename TFixedImage, typename TMovingImage>
virtual const std::vector< ImpactModelConfiguration > & itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetFixedModelsConfiguration ( )
virtual

◆ GetLayersWeight()

template<typename TFixedImage, typename TMovingImage>
virtual std::vector< float > itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetLayersWeight ( ) const
virtual

◆ GetMode()

template<typename TFixedImage, typename TMovingImage>
virtual std::string itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetMode ( ) const
virtual

◆ GetMovingModelsConfiguration()

template<typename TFixedImage, typename TMovingImage>
virtual const std::vector< ImpactModelConfiguration > & itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetMovingModelsConfiguration ( )
virtual

◆ GetPCA()

template<typename TFixedImage, typename TMovingImage>
virtual std::vector< unsigned int > itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetPCA ( ) const
virtual

◆ GetSeed()

template<typename TFixedImage, typename TMovingImage>
virtual unsigned int itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetSeed ( ) const
virtual

◆ GetSubsetFeatures()

template<typename TFixedImage, typename TMovingImage>
virtual std::vector< unsigned int > itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetSubsetFeatures ( ) const
virtual

◆ GetSubsetOfFeatures()

template<typename TFixedImage, typename TMovingImage>
std::vector< unsigned int > itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetSubsetOfFeatures ( const std::vector< unsigned int > & featuresIndex,
std::mt19937 & randomGenerator,
int n ) const
private

Retrieves a subset of features based on the provided indices.

This method selects a subset of features from the full set of features, based on the indices provided in featuresIndex. It uses a random number generator for sampling a subset of the features, ensuring that the selected features are randomly chosen.

Parameters
featuresIndexA vector of indices representing the features to be considered.
randomGeneratorA random number generator used for sampling the features.
nThe number of features to select from the provided list of indices.
Returns
A vector containing the indices of the randomly selected features.

◆ GetUseMixedPrecision()

template<typename TFixedImage, typename TMovingImage>
virtual bool itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetUseMixedPrecision ( ) const
virtual

◆ GetValue()

template<typename TFixedImage, typename TMovingImage>
MeasureType itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetValue ( const ParametersType & parameters) const
override

Compute the similarity value (loss) for a given transformation parameter set. This is the main entry point for single-valued optimizers and is multi-threaded internally. It aggregates the contribution from all threads.

◆ GetValueAndDerivative()

template<typename TFixedImage, typename TMovingImage>
void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetValueAndDerivative ( const ParametersType & parameters,
MeasureType & value,
DerivativeType & derivative ) const
override

Compute both the similarity value and its gradient in a multi-threaded context. This is the main function called by optimizers requiring both value and derivative, and it supports full parallel execution.

◆ GetValueAndDerivativeSingleThreaded()

template<typename TFixedImage, typename TMovingImage>
void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetValueAndDerivativeSingleThreaded ( const ParametersType & parameters,
MeasureType & value,
DerivativeType & derivative ) const

Compute both the similarity value and its gradient in a single-threaded context.

◆ GetValueSingleThreaded()

template<typename TFixedImage, typename TMovingImage>
virtual MeasureType itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetValueSingleThreaded ( const ParametersType & parameters) const
virtual

Compute the similarity value (loss) for a given transformation parameter set. This method is intended for use with single-valued optimizers in a single-threaded context. It is typically used in testing or debugging scenarios.

◆ GetWriteFeatureMaps()

template<typename TFixedImage, typename TMovingImage>
virtual bool itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::GetWriteFeatureMaps ( ) const
virtual

◆ Initialize()

template<typename TFixedImage, typename TMovingImage>
void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::Initialize ( )
override

Initializes the metric and loads models, interpolators, and feature map settings. Called before the optimization loop starts. Ensures all configuration dependencies are resolved.

◆ InitializeThreadingParameters()

template<typename TFixedImage, typename TMovingImage>
void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::InitializeThreadingParameters ( ) const
overrideprotectedvirtual

Initializes per-thread loss structures and ensures thread safety for parallel execution. Overrides superclass method because the metric uses its own loss aggregation system.

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

◆ ITK_DISALLOW_COPY_AND_MOVE()

template<typename TFixedImage, typename TMovingImage>
itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::ITK_DISALLOW_COPY_AND_MOVE ( ImpactImageToImageMetric< TFixedImage, TMovingImage > )

◆ itkAlignedTypedef()

template<typename TFixedImage, typename TMovingImage>
itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::itkAlignedTypedef ( ITK_CACHE_LINE_ALIGNMENT ,
PaddedLossPerThreadStruct ,
AlignedLossPerThreadStruct  )
private

◆ itkPadStruct()

template<typename TFixedImage, typename TMovingImage>
itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::itkPadStruct ( ITK_CACHE_LINE_ALIGNMENT ,
LossPerThreadStruct ,
PaddedLossPerThreadStruct  )
private

Thread-safe wrapper for per-thread loss computation (padded to avoid false sharing).

◆ itkStaticConstMacro() [1/2]

template<typename TFixedImage, typename TMovingImage>
itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::itkStaticConstMacro ( FixedImageDimension ,
unsigned int ,
FixedImageType::ImageDimension  )

The fixed image dimension.

◆ itkStaticConstMacro() [2/2]

template<typename TFixedImage, typename TMovingImage>
itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::itkStaticConstMacro ( MovingImageDimension ,
unsigned int ,
MovingImageType::ImageDimension  )

The moving image dimension.

◆ New()

template<typename TFixedImage, typename TMovingImage>
Pointer itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::New ( )
static

Method for creation through the object factory.

◆ SampleCheck() [1/2]

template<typename TFixedImage, typename TMovingImage>
bool itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SampleCheck ( const FixedImagePointType & fixedImageCenterCoordinate) const
protected

Checks if the fixed image point lies within valid bounds for sampling.

This function checks if a fixed image point is within the valid bounds for sampling. It is used to validate individual points before they are processed in the similarity metric.

Parameters
fixedImageCenterCoordinateThe coordinate of the fixed image point to validate.
Returns
true if the fixed image point is within valid bounds for sampling, false otherwise.

◆ SampleCheck() [2/2]

template<typename TFixedImage, typename TMovingImage>
bool itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SampleCheck ( const FixedImagePointType & fixedImageCenterCoordinate,
const std::vector< std::vector< float > > & patchIndex ) const
protected

Checks if a patch centered at the given fixed image point is valid for sampling.

This function verifies that a patch, defined by the patchIndex, centered at the given fixed image point, is valid for sampling. It ensures that all transformed points of the patch remain inside the domain of the moving image, or within the moving mask if one is defined. It is used to validate patches before they are processed in the similarity metric.

Parameters
fixedImageCenterCoordinateThe center coordinate of the patch in the fixed image.
patchIndexA layout of the patch, defining the indices of the region to sample around the fixed image point.
Returns
true if the patch is valid for sampling, meaning all transformed points are inside the moving image domain or mask; false otherwise.

◆ SetCurrentLevel()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetCurrentLevel ( unsigned int _arg)
virtual

Set/Get the current resolution level

◆ SetDevice()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetDevice ( torch::Device _arg)
virtual

Set/Get the device on which all model inference and tensor operations are performed. Example: torch::Device(torch::kCUDA, 0) for GPU 0.

◆ SetDistance()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetDistance ( std::vector< std::string > _arg)
virtual

Set/Get the type of loss function used for each layer (e.g., "l1", "cosine", "ncc"). Supports heterogeneous losses across layers to adapt to the nature of each feature representation.

◆ SetFeatureMapsPath()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetFeatureMapsPath ( std::string _arg)
virtual

Set/Get the directory path where feature maps will be written if WriteFeatureMaps is true. The path will be created if it does not exist.

◆ SetFeaturesMapUpdateInterval()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetFeaturesMapUpdateInterval ( int _arg)
virtual

Set/Get how often (in number of optimizer iterations) the feature maps should be updated. A value of 0 disables updates (useful in static mode). Positive values enable periodic refreshes.

◆ SetFixedModelsConfiguration()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetFixedModelsConfiguration ( std::vector< ImpactModelConfiguration > _arg)
virtual

Set/Get the list of TorchScript model configurations used to extract features from the fixed image. Each model can target a different resolution, architecture, or semantic level.

◆ SetLayersWeight()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetLayersWeight ( std::vector< float > _arg)
virtual

Set/Get the weights applied to each layer's loss contribution. Useful for balancing the influence of layers with different semantic granularity.

◆ SetMode()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetMode ( std::string _arg)
virtual

Set/Get the mode of operation: "Jacobian", "Static", or "Dynamic".

  • "Jacobian": online patch extraction with gradient backpropagation.
  • "Static": precomputed full feature maps.

◆ SetMovingModelsConfiguration()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetMovingModelsConfiguration ( std::vector< ImpactModelConfiguration > _arg)
virtual

Set/Get the list of TorchScript model configurations used to extract features from the moving image. Allows using different models for fixed and moving images to support asymmetric or multimodal setups.

◆ SetPCA()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetPCA ( std::vector< unsigned int > _arg)
virtual

Set/Get the number of principal components to keep after applying PCA to the feature maps. Set to 0 to disable PCA. Reduces dimensionality and improve runtime.

◆ SetSeed()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetSeed ( unsigned int _arg)
virtual

Set/Get the manual seed

◆ SetSubsetFeatures()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetSubsetFeatures ( std::vector< unsigned int > _arg)
virtual

Set/Get the subset of feature indices to be used in the loss computation. This allows dimensionality reduction or focusing on the most informative channels.

◆ SetUseMixedPrecision()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetUseMixedPrecision ( bool _arg)
virtual

Set/Get whether mixed precision (float16/float32) should be used during model inference.

◆ SetWriteFeatureMaps()

template<typename TFixedImage, typename TMovingImage>
virtual void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::SetWriteFeatureMaps ( bool _arg)
virtual

Set/Get whether the extracted feature maps should be written to disk (for inspection or debugging). Useful for visualizing the intermediate representations used by the metric.

◆ ThreadedGetValue()

template<typename TFixedImage, typename TMovingImage>
void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::ThreadedGetValue ( ThreadIdType threadID) const
overrideprotectedvirtual

Computes the similarity value contribution for a given thread.

This method is called in parallel across multiple threads. Each thread calculates and accumulates its partial loss contribution to the overall similarity value.

Parameters
threadIDThe unique identifier of the thread processing the similarity calculation.

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

◆ ThreadedGetValueAndDerivative()

template<typename TFixedImage, typename TMovingImage>
void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::ThreadedGetValueAndDerivative ( ThreadIdType threadID) const
overrideprotectedvirtual

Computes both the similarity value and its gradient for a given thread.

Each thread computes the semantic loss and its gradient with respect to the transformation parameters.

Parameters
threadIDThe unique identifier of the thread performing the calculation.

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

◆ UpdateFixedFeaturesMaps()

template<typename TFixedImage, typename TMovingImage>
void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::UpdateFixedFeaturesMaps ( )
protected

Updates the fixed feature maps in static mode.

This method re-extracts deep features from the images using the TorchScript models when transitioning to a new pyramid level or after a specified feature update interval. This ensures that the feature maps are kept up-to-date for the registration process.

◆ UpdateMovingFeaturesMaps()

template<typename TFixedImage, typename TMovingImage>
void itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::UpdateMovingFeaturesMaps ( )
protected

Updates the moving feature maps in static mode.

This method performs the same operation as UpdateFixedFeaturesMaps(), but it applies to the moving image. It re-extracts deep features from the moving image using the TorchScript models, ensuring that the feature maps are updated when transitioning to a new pyramid level or after a feature update interval.

Member Data Documentation

◆ m_CurrentLevel

template<typename TFixedImage, typename TMovingImage>
unsigned int itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_CurrentLevel
private

Definition at line 659 of file itkImpactImageToImageMetric.h.

◆ m_Device

template<typename TFixedImage, typename TMovingImage>
torch::Device itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_Device = torch::Device(torch::kCPU)
private

Definition at line 657 of file itkImpactImageToImageMetric.h.

◆ m_Distance

template<typename TFixedImage, typename TMovingImage>
std::vector<std::string> itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_Distance
private

Definition at line 652 of file itkImpactImageToImageMetric.h.

◆ m_FeatureMapsPath

template<typename TFixedImage, typename TMovingImage>
std::string itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_FeatureMapsPath
private

Definition at line 656 of file itkImpactImageToImageMetric.h.

◆ m_FeaturesIndexes

template<typename TFixedImage, typename TMovingImage>
std::vector<std::vector<unsigned int> > itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_FeaturesIndexes
private

Definition at line 667 of file itkImpactImageToImageMetric.h.

◆ m_FeaturesMapUpdateInterval

template<typename TFixedImage, typename TMovingImage>
int itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_FeaturesMapUpdateInterval
private

Definition at line 653 of file itkImpactImageToImageMetric.h.

◆ m_FixedFeaturesMaps

template<typename TFixedImage, typename TMovingImage>
std::vector<FeaturesMaps> itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_FixedFeaturesMaps
private

Definition at line 663 of file itkImpactImageToImageMetric.h.

◆ m_FixedInterpolator

template<typename TFixedImage, typename TMovingImage>
InterpolatorPointer itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_FixedInterpolator
private
Initial value:
= [this] {
const auto interpolator = FixedInterpolatorType::New();
interpolator->SetSplineOrder(3);
return interpolator;
}()

Interpolator for fixed image intensity values, set once at initialization. Uses 3rd-order B-spline interpolation.

Definition at line 674 of file itkImpactImageToImageMetric.h.

◆ m_FixedModelsConfiguration

template<typename TFixedImage, typename TMovingImage>
std::vector<ImpactModelConfiguration> itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_FixedModelsConfiguration
private

TorchScript model configurations for fixed and moving image feature extraction.

Definition at line 646 of file itkImpactImageToImageMetric.h.

◆ m_LayersWeight

template<typename TFixedImage, typename TMovingImage>
std::vector<float> itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_LayersWeight
private

Definition at line 651 of file itkImpactImageToImageMetric.h.

◆ m_LossThreadStruct

template<typename TFixedImage, typename TMovingImage>
std::unique_ptr<AlignedLossPerThreadStruct[]> itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_LossThreadStruct { nullptr }
mutableprivate

Per-thread loss structures, dynamically allocated during initialization.

Definition at line 702 of file itkImpactImageToImageMetric.h.

◆ m_LossThreadStructSize

template<typename TFixedImage, typename TMovingImage>
int itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_LossThreadStructSize = 0
mutableprivate

Definition at line 704 of file itkImpactImageToImageMetric.h.

◆ m_Mode

template<typename TFixedImage, typename TMovingImage>
std::string itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_Mode
private

Definition at line 654 of file itkImpactImageToImageMetric.h.

◆ m_MovingFeaturesMaps

template<typename TFixedImage, typename TMovingImage>
std::vector<FeaturesMaps> itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_MovingFeaturesMaps
private

Definition at line 664 of file itkImpactImageToImageMetric.h.

◆ m_MovingModelsConfiguration

template<typename TFixedImage, typename TMovingImage>
std::vector<ImpactModelConfiguration> itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_MovingModelsConfiguration
private

Definition at line 647 of file itkImpactImageToImageMetric.h.

◆ m_PCA

template<typename TFixedImage, typename TMovingImage>
std::vector<unsigned int> itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_PCA
private

Definition at line 650 of file itkImpactImageToImageMetric.h.

◆ m_PrincipalComponents

template<typename TFixedImage, typename TMovingImage>
std::vector<torch::Tensor> itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_PrincipalComponents
private

Definition at line 665 of file itkImpactImageToImageMetric.h.

◆ m_Seed

template<typename TFixedImage, typename TMovingImage>
unsigned int itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_Seed
private

Definition at line 660 of file itkImpactImageToImageMetric.h.

◆ m_SubsetFeatures

template<typename TFixedImage, typename TMovingImage>
std::vector<unsigned int> itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_SubsetFeatures
private

Definition at line 649 of file itkImpactImageToImageMetric.h.

◆ m_UseMixedPrecision

template<typename TFixedImage, typename TMovingImage>
bool itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_UseMixedPrecision
private

Definition at line 658 of file itkImpactImageToImageMetric.h.

◆ m_WriteFeatureMaps

template<typename TFixedImage, typename TMovingImage>
bool itk::ImpactImageToImageMetric< TFixedImage, TMovingImage >::m_WriteFeatureMaps
private

Definition at line 655 of file itkImpactImageToImageMetric.h.



Generated on 26-02-2026 for elastix by doxygen 1.16.1 (669aeeefca743c148e2d935b3d3c69535c7491e6) elastix logo