SimpleITK  
Public Types | Public Member Functions | Private Types | Private Member Functions | Private Attributes | Friends | List of all members
itk::simple::DisplacementFieldJacobianDeterminantFilter Class Reference

Computes a scalar image from a vector image (e.g., deformation field) input, where each output scalar at each pixel is the Jacobian determinant of the vector field at that location. This calculation is correct in the case where the vector image is a "displacement" from the current location. The computation for the jacobian determinant is: det[ dT/dx ] = det[ I + du/dx ]. More...

#include <sitkDisplacementFieldJacobianDeterminantFilter.h>

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

Public Types

using PixelIDTypeList = RealVectorPixelIDTypeList
 
using Self = DisplacementFieldJacobianDeterminantFilter
 
- Public Types inherited from itk::simple::ImageFilter
using Self = ImageFilter
 
- Public Types inherited from itk::simple::ProcessObject
using Self = ProcessObject
 

Public Member Functions

 DisplacementFieldJacobianDeterminantFilter ()
 
Image Execute (const Image &image1)
 
std::vector< double > GetDerivativeWeights () const
 
std::string GetName () const
 
bool GetUseImageSpacing () const
 
SelfSetDerivativeWeights (std::vector< double > DerivativeWeights)
 
SelfSetUseImageSpacing (bool UseImageSpacing)
 
std::string ToString () const
 
SelfUseImageSpacingOff ()
 
SelfUseImageSpacingOn ()
 
virtual ~DisplacementFieldJacobianDeterminantFilter ()
 
- Public Member Functions inherited from itk::simple::ImageFilter
 ImageFilter ()
 
virtual ~ImageFilter ()=0
 
- 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
 

Private Types

using MemberFunctionType = Image(Self::*)(const Image &image1)
 

Private Member Functions

template<class TImageType >
Image ExecuteInternal (const Image &image1)
 

Private Attributes

std::vector< double > m_DerivativeWeights {std::vector<double>()}
 
std::unique_ptr< detail::MemberFunctionFactory< MemberFunctionType > > m_MemberFactory
 
bool m_UseImageSpacing {true}
 

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...
 
- Protected Member Functions inherited from itk::simple::ImageFilter
void CheckImageMatchingDimension (const Image &image1, const Image &image2, const std::string &image2Name)
 
void CheckImageMatchingPixelType (const Image &image1, const Image &image2, const std::string &image2Name)
 
void CheckImageMatchingSize (const Image &image1, const Image &image2, const std::string &image2Name)
 
- Protected Member Functions inherited from itk::simple::ProcessObject
virtual unsigned long AddITKObserver (const itk::EventObject &, itk::Command *)
 
virtual itk::ProcessObjectGetActiveProcess ()
 
virtual void OnActiveProcessDelete ()
 
virtual void onCommandDelete (const itk::simple::Command *cmd) noexcept
 
virtual void PreUpdate (itk::ProcessObject *p)
 
virtual void RemoveITKObserver (EventCommand &e)
 
- Protected Member Functions inherited from itk::simple::NonCopyable
 NonCopyable ()=default
 
 NonCopyable (const NonCopyable &)=delete
 
NonCopyableoperator= (const NonCopyable &)=delete
 
- Static Protected Member Functions inherited from itk::simple::ImageFilter
template<class TImageType >
static void FixNonZeroIndex (TImageType *img)
 
- 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<unsigned int VImageDimension, unsigned int VLength, template< unsigned int > class TVector>
static Image CastITKToImage (itk::Image< TVector< 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

Computes a scalar image from a vector image (e.g., deformation field) input, where each output scalar at each pixel is the Jacobian determinant of the vector field at that location. This calculation is correct in the case where the vector image is a "displacement" from the current location. The computation for the jacobian determinant is: det[ dT/dx ] = det[ I + du/dx ].

Overview
This filter is based on itkVectorGradientMagnitudeImageFilter and supports the m_DerivativeWeights weights for partial derivatives.

Note that the determinant of a zero vector field is also zero, whereas the Jacobian determinant of the corresponding identity warp transformation is 1.0. In order to compute the effective deformation Jacobian determinant 1.0 must be added to the diagonal elements of Jacobian prior to taking the derivative. i.e. det([ (1.0+dx/dx) dx/dy dx/dz ; dy/dx (1.0+dy/dy) dy/dz; dz/dx dz/dy (1.0+dz/dz) ])

Template Parameters (Input and Output)
This filter has one required template parameter which defines the input image type. The pixel type of the input image is assumed to be a vector (e.g., itk::Vector , itk::RGBPixel , itk::FixedArray ). The scalar type of the vector components must be castable to floating point. Instantiating with an image of RGBPixel<unsigned short>, for example, is allowed, but the filter will convert it to an image of Vector<float,3> for processing.

The second template parameter, TRealType, can be optionally specified to define the scalar numerical type used in calculations. This is the component type of the output image, which will be of itk::Vector<TRealType, N>, where N is the number of channels in the multiple component input image. The default type of TRealType is float. For extra precision, you may safely change this parameter to double.

The third template parameter is the output image type. The third parameter will be automatically constructed from the first and second parameters, so it is not necessary (or advisable) to set this parameter explicitly. Given an M-channel input image with dimensionality N, and a numerical type specified as TRealType, the output image will be of type itk::Image<TRealType, N>.

Filter Parameters
The method UseImageSpacingOn will cause derivatives in the image to be scaled (inversely) with the pixel size of the input image, effectively taking derivatives in world coordinates (versus isotropic image space). UseImageSpacingOff turns this functionality off. Default is UseImageSpacingOn. The parameter UseImageSpacing can be set directly with the method SetUseImageSpacing(bool) .

Weights can be applied to the derivatives directly using the SetDerivativeWeights method. Note that if UseImageSpacing is set to TRUE (ON), then these weights will be overridden by weights derived from the image spacing when the filter is updated. The argument to this method is a C array of TRealValue type.

Constraints
We use vnl_det for determinant computation, which only supports square matrices. So the vector dimension of the input image values must be equal to the image dimensions, which is trivially true for a deformation field that maps an n-dimensional space onto itself.

Currently, dimensions up to and including 4 are supported. This limitation comes from the presence of vnl_det() functions for matrices of dimension up to 4x4.

The template parameter TRealType must be floating point (float or double) or a user-defined "real" numerical type with arithmetic operations defined sufficient to compute derivatives.

See also
Image
Neighborhood
NeighborhoodOperator
NeighborhoodIterator
Note
This class was adapted by
Author
Hans J. Johnson, The University of Iowa from code provided by
Tom Vercauteren, INRIA & Mauna Kea Technologies
Torsten Rohlfing, Neuroscience Program, SRI International.
See also
itk::simple::DisplacementFieldJacobianDeterminantFilter for the procedural interface
itk::DisplacementFieldJacobianDeterminantFilter for the Doxygen on the original ITK class.

Definition at line 90 of file sitkDisplacementFieldJacobianDeterminantFilter.h.

Member Typedef Documentation

◆ MemberFunctionType

Setup for member function dispatching

Definition at line 144 of file sitkDisplacementFieldJacobianDeterminantFilter.h.

◆ PixelIDTypeList

Define the pixels types supported by this filter

Definition at line 102 of file sitkDisplacementFieldJacobianDeterminantFilter.h.

◆ Self

Constructor & Destructor Documentation

◆ ~DisplacementFieldJacobianDeterminantFilter()

virtual itk::simple::DisplacementFieldJacobianDeterminantFilter::~DisplacementFieldJacobianDeterminantFilter ( )
virtual

Destructor

◆ DisplacementFieldJacobianDeterminantFilter()

itk::simple::DisplacementFieldJacobianDeterminantFilter::DisplacementFieldJacobianDeterminantFilter ( )

Default Constructor that takes no arguments and initializes default parameters

Member Function Documentation

◆ Execute()

Image itk::simple::DisplacementFieldJacobianDeterminantFilter::Execute ( const Image image1)

Execute the filter on the input image

◆ ExecuteInternal()

template<class TImageType >
Image itk::simple::DisplacementFieldJacobianDeterminantFilter::ExecuteInternal ( const Image image1)
private

◆ GetDerivativeWeights()

std::vector<double> itk::simple::DisplacementFieldJacobianDeterminantFilter::GetDerivativeWeights ( ) const
inline

Directly Set/Get the array of weights used in the gradient calculations. Note that calling UseImageSpacingOn will clobber these values.

Definition at line 127 of file sitkDisplacementFieldJacobianDeterminantFilter.h.

◆ GetName()

std::string itk::simple::DisplacementFieldJacobianDeterminantFilter::GetName ( ) const
inlinevirtual

Name of this class

Implements itk::simple::ProcessObject.

Definition at line 130 of file sitkDisplacementFieldJacobianDeterminantFilter.h.

◆ GetUseImageSpacing()

bool itk::simple::DisplacementFieldJacobianDeterminantFilter::GetUseImageSpacing ( ) const
inline

Set/Get whether or not the filter will use the spacing of the input image (1/spacing) in the calculation of the Jacobian determinant. Use On to compute the Jacobian determinant in the space in which the data was acquired; use Off to reset the derivative weights, ignore the image spacing, and to compute the Jacobian determinant in the image space. Default is On.

Definition at line 117 of file sitkDisplacementFieldJacobianDeterminantFilter.h.

◆ SetDerivativeWeights()

Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::SetDerivativeWeights ( std::vector< double >  DerivativeWeights)
inline

Directly Set/Get the array of weights used in the gradient calculations. Note that calling UseImageSpacingOn will clobber these values.

Definition at line 122 of file sitkDisplacementFieldJacobianDeterminantFilter.h.

◆ SetUseImageSpacing()

Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::SetUseImageSpacing ( bool  UseImageSpacing)
inline

Set/Get whether or not the filter will use the spacing of the input image (1/spacing) in the calculation of the Jacobian determinant. Use On to compute the Jacobian determinant in the space in which the data was acquired; use Off to reset the derivative weights, ignore the image spacing, and to compute the Jacobian determinant in the image space. Default is On.

Definition at line 108 of file sitkDisplacementFieldJacobianDeterminantFilter.h.

◆ ToString()

std::string itk::simple::DisplacementFieldJacobianDeterminantFilter::ToString ( ) const
virtual

Print ourselves out

Reimplemented from itk::simple::ProcessObject.

◆ UseImageSpacingOff()

Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::UseImageSpacingOff ( )
inline

◆ UseImageSpacingOn()

Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::UseImageSpacingOn ( )
inline

Set the value of UseImageSpacing to true or false respectfully.

Definition at line 111 of file sitkDisplacementFieldJacobianDeterminantFilter.h.

Friends And Related Function Documentation

◆ detail::MemberFunctionAddressor< MemberFunctionType >

Member Data Documentation

◆ m_DerivativeWeights

std::vector<double> itk::simple::DisplacementFieldJacobianDeterminantFilter::m_DerivativeWeights {std::vector<double>()}
private

◆ m_MemberFactory

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

◆ m_UseImageSpacing

bool itk::simple::DisplacementFieldJacobianDeterminantFilter::m_UseImageSpacing {true}
private

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