Package org.itk.simple
Class DisplacementFieldJacobianDeterminantFilter
java.lang.Object
org.itk.simple.ProcessObject
org.itk.simple.ImageFilter
org.itk.simple.DisplacementFieldJacobianDeterminantFilter
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:
Image
Neighborhood
NeighborhoodOperator
NeighborhoodIterator
This class was adapted by
Hans J. Johnson, The University of Iowa from code provided by
Tom Vercauteren, INRIA & Mauna Kea Technologies
Torsten Rohlfing, Neuroscience Program, SRI International.
See:
itk::simple::DisplacementFieldJacobianDeterminantFilter for the procedural interface
itk::DisplacementFieldJacobianDeterminantFilter for the Doxygen on the original ITK class.
C++ includes: sitkDisplacementFieldJacobianDeterminantFilter.h
-
Field Summary
Fields inherited from class org.itk.simple.ProcessObject
swigCMemOwn
-
Constructor Summary
ConstructorsModifierConstructorDescriptionitk::simple::DisplacementFieldJacobianDeterminantFilter::DisplacementFieldJacobianDeterminantFilter() Default Constructor that takes no arguments and initializes default parametersprotected
DisplacementFieldJacobianDeterminantFilter
(long cPtr, boolean cMemoryOwn) -
Method Summary
Modifier and TypeMethodDescriptionvoid
delete()
virtual itk::simple::DisplacementFieldJacobianDeterminantFilter::~DisplacementFieldJacobianDeterminantFilter() DestructorImage itk::simple::DisplacementFieldJacobianDeterminantFilter::Execute(const Image &image1) Execute the filter on the input imageprotected void
finalize()
protected static long
std::vector<double> itk::simple::DisplacementFieldJacobianDeterminantFilter::GetDerivativeWeights() const Directly Set/Get the array of weights used in the gradient calculations.getName()
std::string itk::simple::DisplacementFieldJacobianDeterminantFilter::GetName() const Name of this classboolean
bool itk::simple::DisplacementFieldJacobianDeterminantFilter::GetUseImageSpacing() const Set/Get whether or not the filter will use the spacing of the input image (1/spacing) in the calculation of the Jacobian determinant.void
setDerivativeWeights
(VectorDouble DerivativeWeights) Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::SetDerivativeWeights(std::vector< double > DerivativeWeights) Directly Set/Get the array of weights used in the gradient calculations.void
setUseImageSpacing
(boolean UseImageSpacing) Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::SetUseImageSpacing(bool UseImageSpacing) Set/Get whether or not the filter will use the spacing of the input image (1/spacing) in the calculation of the Jacobian determinant.protected static long
toString()
std::string itk::simple::DisplacementFieldJacobianDeterminantFilter::ToString() const Print ourselves outvoid
Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::UseImageSpacingOff()void
Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::UseImageSpacingOn() Set the value of UseImageSpacing to true or false respectfully.Methods inherited from class org.itk.simple.ImageFilter
getCPtr, swigRelease
Methods inherited from class org.itk.simple.ProcessObject
abort, addCommand, debugOff, debugOn, getCPtr, getDebug, getGlobalDefaultCoordinateTolerance, getGlobalDefaultDebug, getGlobalDefaultDirectionTolerance, getGlobalDefaultNumberOfThreads, getGlobalDefaultThreader, getGlobalWarningDisplay, getNumberOfThreads, getNumberOfWorkUnits, getProgress, globalDefaultDebugOff, globalDefaultDebugOn, globalWarningDisplayOff, globalWarningDisplayOn, hasCommand, removeAllCommands, setDebug, setGlobalDefaultCoordinateTolerance, setGlobalDefaultDebug, setGlobalDefaultDirectionTolerance, setGlobalDefaultNumberOfThreads, setGlobalDefaultThreader, setGlobalWarningDisplay, setNumberOfThreads, setNumberOfWorkUnits, swigRelease
-
Constructor Details
-
DisplacementFieldJacobianDeterminantFilter
protected DisplacementFieldJacobianDeterminantFilter(long cPtr, boolean cMemoryOwn) -
DisplacementFieldJacobianDeterminantFilter
public DisplacementFieldJacobianDeterminantFilter()itk::simple::DisplacementFieldJacobianDeterminantFilter::DisplacementFieldJacobianDeterminantFilter() Default Constructor that takes no arguments and initializes default parameters
-
-
Method Details
-
getCPtr
-
swigRelease
-
finalize
protected void finalize()- Overrides:
finalize
in classImageFilter
-
delete
public void delete()virtual itk::simple::DisplacementFieldJacobianDeterminantFilter::~DisplacementFieldJacobianDeterminantFilter() Destructor- Overrides:
delete
in classImageFilter
-
setUseImageSpacing
public void setUseImageSpacing(boolean UseImageSpacing) Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::SetUseImageSpacing(bool UseImageSpacing) 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. -
useImageSpacingOn
public void useImageSpacingOn()Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::UseImageSpacingOn() Set the value of UseImageSpacing to true or false respectfully. -
useImageSpacingOff
public void useImageSpacingOff()Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::UseImageSpacingOff() -
getUseImageSpacing
public boolean getUseImageSpacing()bool itk::simple::DisplacementFieldJacobianDeterminantFilter::GetUseImageSpacing() const 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. -
setDerivativeWeights
Self& itk::simple::DisplacementFieldJacobianDeterminantFilter::SetDerivativeWeights(std::vector< double > DerivativeWeights) Directly Set/Get the array of weights used in the gradient calculations. Note that calling UseImageSpacingOn will clobber these values. -
getDerivativeWeights
std::vector<double> itk::simple::DisplacementFieldJacobianDeterminantFilter::GetDerivativeWeights() const Directly Set/Get the array of weights used in the gradient calculations. Note that calling UseImageSpacingOn will clobber these values. -
getName
std::string itk::simple::DisplacementFieldJacobianDeterminantFilter::GetName() const Name of this class- Overrides:
getName
in classProcessObject
-
toString
std::string itk::simple::DisplacementFieldJacobianDeterminantFilter::ToString() const Print ourselves out- Overrides:
toString
in classProcessObject
-
execute
Image itk::simple::DisplacementFieldJacobianDeterminantFilter::Execute(const Image &image1) Execute the filter on the input image
-