SimpleITK
1.0.1
|
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>
Public Types | |
typedef RealVectorPixelIDTypeList | PixelIDTypeList |
typedef DisplacementFieldJacobianDeterminantFilter | Self |
Public Types inherited from itk::simple::ImageFilter< 1 > | |
typedef ImageFilter | Self |
Public Types inherited from itk::simple::ProcessObject | |
typedef ProcessObject | Self |
Public Member Functions | |
DisplacementFieldJacobianDeterminantFilter () | |
Image | Execute (const Image &image1) |
Image | Execute (const Image &image1, bool useImageSpacing, const std::vector< double > &derivativeWeights) |
std::vector< double > | GetDerivativeWeights () const |
std::string | GetName () const |
bool | GetUseImageSpacing () const |
Self & | SetDerivativeWeights (const std::vector< double > &DerivativeWeights) |
Self & | SetUseImageSpacing (bool UseImageSpacing) |
std::string | ToString () const |
Self & | UseImageSpacingOff () |
Self & | UseImageSpacingOn () |
~DisplacementFieldJacobianDeterminantFilter () | |
Public Member Functions inherited from itk::simple::ImageFilter< 1 > | |
ImageFilter () | |
virtual | ~ImageFilter ()=0 |
Public Member Functions inherited from itk::simple::ProcessObject | |
virtual void | Abort () |
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 |
Private Types | |
typedef Image(Self::* | MemberFunctionType) (const Image &image1) |
Private Member Functions | |
template<class TImageType > | |
Image | ExecuteInternal (const Image &image1) |
Private Attributes | |
std::vector< double > | m_DerivativeWeights |
nsstd::auto_ptr< detail::MemberFunctionFactory< MemberFunctionType > > | m_MemberFactory |
bool | m_UseImageSpacing |
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 void | SetGlobalDefaultNumberOfThreads (unsigned int n) |
static unsigned int | GetGlobalDefaultNumberOfThreads () |
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... | |
Protected Member Functions inherited from itk::simple::ProcessObject | |
virtual unsigned long | AddITKObserver (const itk::EventObject &, itk::Command *) |
virtual itk::ProcessObject * | GetActiveProcess () |
virtual void | OnActiveProcessDelete () |
virtual void | onCommandDelete (const itk::simple::Command *cmd) throw () |
virtual void | PreUpdate (itk::ProcessObject *p) |
virtual void | RemoveITKObserver (EventCommand &e) |
Protected Member Functions inherited from itk::simple::NonCopyable | |
NonCopyable () | |
Static Protected Member Functions inherited from itk::simple::ImageFilter< 1 > | |
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 TImageType > | |
static Image | CastITKToImage (TImageType *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) |
static const itk::EventObject & | GetITKEventObject (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) |
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 ].
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) ])
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>.
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.
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.
Tom Vercauteren, INRIA & Mauna Kea Technologies
Torsten Rohlfing, Neuroscience Program, SRI International.
Definition at line 80 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
|
private |
Setup for member function dispatching
Definition at line 137 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
typedef RealVectorPixelIDTypeList itk::simple::DisplacementFieldJacobianDeterminantFilter::PixelIDTypeList |
Define the pixels types supported by this filter
Definition at line 92 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
typedef DisplacementFieldJacobianDeterminantFilter itk::simple::DisplacementFieldJacobianDeterminantFilter::Self |
Definition at line 82 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
itk::simple::DisplacementFieldJacobianDeterminantFilter::DisplacementFieldJacobianDeterminantFilter | ( | ) |
Default Constructor that takes no arguments and initializes default parameters
itk::simple::DisplacementFieldJacobianDeterminantFilter::~DisplacementFieldJacobianDeterminantFilter | ( | ) |
Destructor
Execute the filter on the input image
Image itk::simple::DisplacementFieldJacobianDeterminantFilter::Execute | ( | const Image & | image1, |
bool | useImageSpacing, | ||
const std::vector< double > & | derivativeWeights | ||
) |
Execute the filter on the input image with the given parameters
|
private |
|
inline |
Directly Set/Get the array of weights used in the gradient calculations. Note that calling UseImageSpacingOn will clobber these values.
Definition at line 117 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
|
inlinevirtual |
Name of this class
Implements itk::simple::ProcessObject.
Definition at line 119 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
|
inline |
Definition at line 107 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
|
inline |
Directly Set/Get the array of weights used in the gradient calculations. Note that calling UseImageSpacingOn will clobber these values.
Definition at line 112 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
|
inline |
Set/Get whether or not the filter will use the spacing of the input image in its calculations
Definition at line 99 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
|
virtual |
Print ourselves out
Reimplemented from itk::simple::ProcessObject.
|
inline |
Definition at line 103 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
|
inline |
Set the value of UseImageSpacing to true or false respectfully.
Definition at line 102 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
|
friend |
Definition at line 142 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
|
private |
Definition at line 148 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
|
private |
Definition at line 144 of file sitkDisplacementFieldJacobianDeterminantFilter.h.
|
private |
Definition at line 147 of file sitkDisplacementFieldJacobianDeterminantFilter.h.