ITK  6.0.0
Insight Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage > Class Template Reference

#include <itkDisplacementFieldJacobianDeterminantFilter.h>

Detailed Description

template<typename TInputImage, typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
class itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >

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.

Definition at line 114 of file itkDisplacementFieldJacobianDeterminantFilter.h.

+ Inheritance diagram for itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >:
+ Collaboration diagram for itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >:

Public Types

using ConstNeighborhoodIteratorType = ConstNeighborhoodIterator< RealVectorImageType >
 
using ConstPointer = SmartPointer< const Self >
 
using InputImagePointer = typename InputImageType::Pointer
 
using InputImageType = TInputImage
 
using InputPixelType = typename TInputImage::PixelType
 
using OutputImagePointer = typename OutputImageType::Pointer
 
using OutputImageType = TOutputImage
 
using OutputPixelType = typename TOutputImage::PixelType
 
using Pointer = SmartPointer< Self >
 
using RadiusType = typename ConstNeighborhoodIteratorType::RadiusType
 
using RealType = TRealType
 
using RealVectorImageType = Image< RealVectorType, TInputImage::ImageDimension >
 
using RealVectorType = Vector< TRealType, InputPixelType::Dimension >
 
using Self = DisplacementFieldJacobianDeterminantFilter
 
using Superclass = ImageToImageFilter< TInputImage, TOutputImage >
 
using WeightsType = FixedArray< TRealType, ImageDimension >
 
- Public Types inherited from itk::ImageToImageFilter< TInputImage, TOutputImage >
using ConstPointer = SmartPointer< const Self >
 
using InputImageConstPointer = typename InputImageType::ConstPointer
 
using InputImagePixelType = typename InputImageType::PixelType
 
using InputImagePointer = typename InputImageType::Pointer
 
using InputImageRegionType = typename InputImageType::RegionType
 
using InputImageType = TInputImage
 
using Pointer = SmartPointer< Self >
 
using Self = ImageToImageFilter
 
using Superclass = ImageSource< TOutputImage >
 
- Public Types inherited from itk::ImageSource< TOutputImage >
using ConstPointer = SmartPointer< const Self >
 
using DataObjectIdentifierType = Superclass::DataObjectIdentifierType
 
using DataObjectPointer = DataObject::Pointer
 
using DataObjectPointerArraySizeType = Superclass::DataObjectPointerArraySizeType
 
using OutputImagePixelType = typename OutputImageType::PixelType
 
using OutputImagePointer = typename OutputImageType::Pointer
 
using OutputImageRegionType = typename OutputImageType::RegionType
 
using OutputImageType = TOutputImage
 
using Pointer = SmartPointer< Self >
 
using Self = ImageSource
 
using Superclass = ProcessObject
 
- Public Types inherited from itk::ProcessObject
using ConstPointer = SmartPointer< const Self >
 
using DataObjectIdentifierType = DataObject::DataObjectIdentifierType
 
using DataObjectPointer = DataObject::Pointer
 
using DataObjectPointerArray = std::vector< DataObjectPointer >
 
using DataObjectPointerArraySizeType = DataObjectPointerArray::size_type
 
using MultiThreaderType = MultiThreaderBase
 
using NameArray = std::vector< DataObjectIdentifierType >
 
using Pointer = SmartPointer< Self >
 
using Self = ProcessObject
 
using Superclass = Object
 
- Public Types inherited from itk::Object
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = Object
 
using Superclass = LightObject
 
- Public Types inherited from itk::LightObject
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using Self = LightObject
 

Public Member Functions

void GenerateInputRequestedRegion () override
 
const char * GetNameOfClass () const override
 
void SetUseImageSpacing (bool)
 
virtual bool GetUseImageSpacing () const
 
virtual void UseImageSpacingOn ()
 
void SetDerivativeWeights (const WeightsType &)
 
virtual const WeightsTypeGetDerivativeWeights () const
 
- Public Member Functions inherited from itk::ImageToImageFilter< TInputImage, TOutputImage >
const InputImageTypeGetInput () const
 
const InputImageTypeGetInput (unsigned int idx) const
 
void PopBackInput () override
 
void PopFrontInput () override
 
virtual void PushBackInput (const InputImageType *input)
 
virtual void PushFrontInput (const InputImageType *input)
 
virtual void SetInput (const DataObjectIdentifierType &key, DataObject *input)
 
virtual void SetInput (const InputImageType *input)
 
virtual void SetInput (unsigned int, const TInputImage *image)
 
virtual void SetCoordinateTolerance (double _arg)
 
virtual double GetCoordinateTolerance () const
 
virtual void SetDirectionTolerance (double _arg)
 
virtual double GetDirectionTolerance () const
 
- Public Member Functions inherited from itk::ImageSource< TOutputImage >
OutputImageTypeGetOutput (unsigned int idx)
 
OutputImageTypeGetOutput ()
 
const OutputImageTypeGetOutput () const
 
virtual void GraftOutput (DataObject *graft)
 
virtual void GraftOutput (const DataObjectIdentifierType &key, DataObject *graft)
 
virtual void GraftNthOutput (unsigned int idx, DataObject *graft)
 
ProcessObject::DataObjectPointer MakeOutput (ProcessObject::DataObjectPointerArraySizeType idx) override
 
ProcessObject::DataObjectPointer MakeOutput (const ProcessObject::DataObjectIdentifierType &) override
 
- Public Member Functions inherited from itk::ProcessObject
virtual void AbortGenerateDataOn ()
 
virtual void EnlargeOutputRequestedRegion (DataObject *)
 
virtual const bool & GetAbortGenerateData () const
 
DataObjectPointerArray GetIndexedInputs ()
 
DataObjectPointerArray GetIndexedOutputs ()
 
NameArray GetInputNames () const
 
DataObjectPointerArray GetInputs ()
 
MultiThreaderTypeGetMultiThreader () const
 
DataObjectPointerArraySizeType GetNumberOfIndexedInputs () const
 
DataObjectPointerArraySizeType GetNumberOfIndexedOutputs () const
 
DataObjectPointerArraySizeType GetNumberOfInputs () const
 
DataObjectPointerArraySizeType GetNumberOfOutputs () const
 
virtual DataObjectPointerArraySizeType GetNumberOfValidRequiredInputs () const
 
NameArray GetOutputNames () const
 
DataObjectPointerArray GetOutputs ()
 
virtual float GetProgress () const
 
NameArray GetRequiredInputNames () const
 
bool HasInput (const DataObjectIdentifierType &key) const
 
bool HasOutput (const DataObjectIdentifierType &key) const
 
void IncrementProgress (float increment)
 
virtual void PrepareOutputs ()
 
virtual void PropagateRequestedRegion (DataObject *output)
 
virtual void ResetPipeline ()
 
virtual void SetAbortGenerateData (bool _arg)
 
void SetMultiThreader (MultiThreaderType *threader)
 
virtual void Update ()
 
virtual void UpdateLargestPossibleRegion ()
 
virtual void UpdateOutputData (DataObject *output)
 
virtual void UpdateOutputInformation ()
 
void UpdateProgress (float progress)
 
virtual void SetReleaseDataFlag (bool val)
 
virtual bool GetReleaseDataFlag () const
 
void ReleaseDataFlagOn ()
 
void ReleaseDataFlagOff ()
 
virtual void SetReleaseDataBeforeUpdateFlag (bool _arg)
 
virtual const bool & GetReleaseDataBeforeUpdateFlag () const
 
virtual void ReleaseDataBeforeUpdateFlagOn ()
 
virtual void SetNumberOfWorkUnits (ThreadIdType _arg)
 
virtual const ThreadIdTypeGetNumberOfWorkUnits () const
 
- Public Member Functions inherited from itk::Object
unsigned long AddObserver (const EventObject &event, Command *cmd) const
 
unsigned long AddObserver (const EventObject &event, std::function< void(const EventObject &)> function) const
 
LightObject::Pointer CreateAnother () const override
 
virtual void DebugOff () const
 
virtual void DebugOn () const
 
CommandGetCommand (unsigned long tag)
 
bool GetDebug () const
 
MetaDataDictionaryGetMetaDataDictionary ()
 
const MetaDataDictionaryGetMetaDataDictionary () const
 
virtual ModifiedTimeType GetMTime () const
 
virtual const TimeStampGetTimeStamp () const
 
bool HasObserver (const EventObject &event) const
 
void InvokeEvent (const EventObject &)
 
void InvokeEvent (const EventObject &) const
 
virtual void Modified () const
 
void Register () const override
 
void RemoveAllObservers ()
 
void RemoveObserver (unsigned long tag) const
 
void SetDebug (bool debugFlag) const
 
void SetReferenceCount (int) override
 
void UnRegister () const noexcept override
 
void SetMetaDataDictionary (const MetaDataDictionary &rhs)
 
void SetMetaDataDictionary (MetaDataDictionary &&rrhs)
 
virtual void SetObjectName (std::string _arg)
 
virtual const std::string & GetObjectName () const
 
- Public Member Functions inherited from itk::LightObject
Pointer Clone () const
 
virtual void Delete ()
 
virtual int GetReferenceCount () const
 
void Print (std::ostream &os, Indent indent=0) const
 

Static Public Member Functions

static Pointer New ()
 
- Static Public Member Functions inherited from itk::ImageToImageFilter< TInputImage, TOutputImage >
static double GetGlobalDefaultCoordinateTolerance ()
 
static double GetGlobalDefaultDirectionTolerance ()
 
static void SetGlobalDefaultCoordinateTolerance (double)
 
static void SetGlobalDefaultDirectionTolerance (double)
 
- Static Public Member Functions inherited from itk::Object
static bool GetGlobalWarningDisplay ()
 
static void GlobalWarningDisplayOff ()
 
static void GlobalWarningDisplayOn ()
 
static Pointer New ()
 
static void SetGlobalWarningDisplay (bool val)
 
- Static Public Member Functions inherited from itk::LightObject
static void BreakOnError ()
 
static Pointer New ()
 

Static Public Attributes

static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension
 
static constexpr unsigned int VectorDimension = InputPixelType::Dimension
 
- Static Public Attributes inherited from itk::ImageToImageFilter< TInputImage, TOutputImage >
static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension
 
static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension
 
- Static Public Attributes inherited from itk::ImageSource< TOutputImage >
static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension
 

Protected Types

using ImageBaseType = typename InputImageType::Superclass
 
- Protected Types inherited from itk::ImageToImageFilter< TInputImage, TOutputImage >
using InputToOutputRegionCopierType = ImageToImageFilterDetail::ImageRegionCopier< Self::OutputImageDimension, Self::InputImageDimension >
 
using OutputToInputRegionCopierType = ImageToImageFilterDetail::ImageRegionCopier< Self::InputImageDimension, Self::OutputImageDimension >
 

Protected Member Functions

void BeforeThreadedGenerateData () override
 
 DisplacementFieldJacobianDeterminantFilter ()
 
void DynamicThreadedGenerateData (const OutputImageRegionType &outputRegionForThread) override
 
virtual TRealType EvaluateAtNeighborhood (const ConstNeighborhoodIteratorType &it) const
 
virtual const ImageBaseTypeGetRealValuedInputImage () const
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
 ~DisplacementFieldJacobianDeterminantFilter () override=default
 
virtual const RadiusTypeGetNeighborhoodRadius () const
 
virtual void SetNeighborhoodRadius (RadiusType _arg)
 
- Protected Member Functions inherited from itk::ImageToImageFilter< TInputImage, TOutputImage >
virtual void CallCopyInputRegionToOutputRegion (OutputImageRegionType &destRegion, const InputImageRegionType &srcRegion)
 
virtual void CallCopyOutputRegionToInputRegion (InputImageRegionType &destRegion, const OutputImageRegionType &srcRegion)
 
 ImageToImageFilter ()
 
void VerifyInputInformation () const override
 
 ~ImageToImageFilter () override=default
 
virtual void PushBackInput (const DataObject *input)
 
virtual void PushFrontInput (const DataObject *input)
 
- Protected Member Functions inherited from itk::ImageSource< TOutputImage >
virtual void AfterThreadedGenerateData ()
 
virtual void AllocateOutputs ()
 
void ClassicMultiThread (ThreadFunctionType callbackFunction)
 
void GenerateData () override
 
virtual const ImageRegionSplitterBaseGetImageRegionSplitter () const
 
 ImageSource ()
 
virtual unsigned int SplitRequestedRegion (unsigned int i, unsigned int pieces, OutputImageRegionType &splitRegion)
 
 ~ImageSource () override=default
 
virtual void ThreadedGenerateData (const OutputImageRegionType &region, ThreadIdType threadId)
 
virtual bool GetDynamicMultiThreading () const
 
virtual void SetDynamicMultiThreading (bool _arg)
 
virtual void DynamicMultiThreadingOn ()
 
- Protected Member Functions inherited from itk::ProcessObject
virtual void AddInput (DataObject *input)
 
void AddOptionalInputName (const DataObjectIdentifierType &)
 
void AddOptionalInputName (const DataObjectIdentifierType &, DataObjectPointerArraySizeType idx)
 
virtual void AddOutput (DataObject *output)
 
bool AddRequiredInputName (const DataObjectIdentifierType &)
 
bool AddRequiredInputName (const DataObjectIdentifierType &, DataObjectPointerArraySizeType idx)
 
virtual void CacheInputReleaseDataFlags ()
 
virtual void GenerateOutputInformation ()
 
virtual void GenerateOutputRequestedRegion (DataObject *output)
 
DataObjectGetInput (const DataObjectIdentifierType &key)
 
const DataObjectGetInput (const DataObjectIdentifierType &key) const
 
virtual const DataObjectPointerArraySizeTypeGetNumberOfRequiredInputs () const
 
virtual const DataObjectPointerArraySizeTypeGetNumberOfRequiredOutputs () const
 
bool IsIndexedInputName (const DataObjectIdentifierType &) const
 
bool IsIndexedOutputName (const DataObjectIdentifierType &) const
 
bool IsRequiredInputName (const DataObjectIdentifierType &) const
 
DataObjectPointerArraySizeType MakeIndexFromInputName (const DataObjectIdentifierType &name) const
 
DataObjectPointerArraySizeType MakeIndexFromOutputName (const DataObjectIdentifierType &name) const
 
DataObjectIdentifierType MakeNameFromInputIndex (DataObjectPointerArraySizeType idx) const
 
DataObjectIdentifierType MakeNameFromOutputIndex (DataObjectPointerArraySizeType idx) const
 
 ProcessObject ()
 
virtual void PropagateResetPipeline ()
 
virtual void PushBackInput (const DataObject *input)
 
virtual void PushFrontInput (const DataObject *input)
 
virtual void ReleaseInputs ()
 
virtual void RemoveInput (const DataObjectIdentifierType &key)
 
virtual void RemoveInput (DataObjectPointerArraySizeType)
 
virtual void RemoveOutput (const DataObjectIdentifierType &key)
 
virtual void RemoveOutput (DataObjectPointerArraySizeType idx)
 
bool RemoveRequiredInputName (const DataObjectIdentifierType &)
 
virtual void RestoreInputReleaseDataFlags ()
 
virtual void SetInput (const DataObjectIdentifierType &key, DataObject *input)
 
virtual void SetNthInput (DataObjectPointerArraySizeType idx, DataObject *input)
 
virtual void SetNthOutput (DataObjectPointerArraySizeType idx, DataObject *output)
 
void SetNumberOfIndexedInputs (DataObjectPointerArraySizeType num)
 
void SetNumberOfIndexedOutputs (DataObjectPointerArraySizeType num)
 
virtual void SetNumberOfRequiredInputs (DataObjectPointerArraySizeType)
 
virtual void SetNumberOfRequiredOutputs (DataObjectPointerArraySizeType _arg)
 
virtual void SetOutput (const DataObjectIdentifierType &name, DataObject *output)
 
virtual void SetPrimaryInput (DataObject *object)
 
virtual void SetPrimaryOutput (DataObject *object)
 
void SetRequiredInputNames (const NameArray &)
 
virtual void VerifyPreconditions () const
 
 ~ProcessObject () override
 
DataObjectGetInput (DataObjectPointerArraySizeType idx)
 
const DataObjectGetInput (DataObjectPointerArraySizeType idx) const
 
DataObjectGetPrimaryInput ()
 
const DataObjectGetPrimaryInput () const
 
virtual void SetPrimaryInputName (const DataObjectIdentifierType &key)
 
virtual const char * GetPrimaryInputName () const
 
DataObjectGetOutput (const DataObjectIdentifierType &key)
 
const DataObjectGetOutput (const DataObjectIdentifierType &key) const
 
virtual void SetPrimaryOutputName (const DataObjectIdentifierType &key)
 
virtual const char * GetPrimaryOutputName () const
 
DataObjectGetOutput (DataObjectPointerArraySizeType i)
 
const DataObjectGetOutput (DataObjectPointerArraySizeType i) const
 
DataObjectGetPrimaryOutput ()
 
const DataObjectGetPrimaryOutput () const
 
virtual bool GetThreaderUpdateProgress () const
 
virtual void ThreaderUpdateProgressOn ()
 
virtual void SetThreaderUpdateProgress (bool arg)
 
- Protected Member Functions inherited from itk::Object
 Object ()
 
bool PrintObservers (std::ostream &os, Indent indent) const
 
virtual void SetTimeStamp (const TimeStamp &timeStamp)
 
 ~Object () override
 
- Protected Member Functions inherited from itk::LightObject
virtual LightObject::Pointer InternalClone () const
 
 LightObject ()
 
virtual void PrintHeader (std::ostream &os, Indent indent) const
 
virtual void PrintTrailer (std::ostream &os, Indent indent) const
 
virtual ~LightObject ()
 

Protected Attributes

WeightsType m_DerivativeWeights {}
 
WeightsType m_HalfDerivativeWeights {}
 
- Protected Attributes inherited from itk::ImageSource< TOutputImage >
bool m_DynamicMultiThreading { true }
 
- Protected Attributes inherited from itk::ProcessObject
TimeStamp m_OutputInformationMTime {}
 
bool m_Updating {}
 
- Protected Attributes inherited from itk::LightObject
std::atomic< int > m_ReferenceCount {}
 

Private Attributes

RadiusType m_NeighborhoodRadius {}
 
ImageBaseType::ConstPointer m_RealValuedInputImage {}
 
ThreadIdType m_RequestedNumberOfWorkUnits {}
 
bool m_UseImageSpacing { true }
 

Additional Inherited Members

- Static Protected Member Functions inherited from itk::ImageSource< TOutputImage >
static const ImageRegionSplitterBaseGetGlobalDefaultSplitter ()
 
static ITK_THREAD_RETURN_FUNCTION_CALL_CONVENTION ThreaderCallback (void *arg)
 
- Static Protected Member Functions inherited from itk::ProcessObject
template<typename TSourceObject >
static void MakeRequiredOutputs (TSourceObject &sourceObject, const DataObjectPointerArraySizeType numberOfRequiredOutputs)
 
static constexpr float progressFixedToFloat (uint32_t fixed)
 
static uint32_t progressFloatToFixed (float f)
 

Member Typedef Documentation

◆ ConstNeighborhoodIteratorType

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::ConstNeighborhoodIteratorType = ConstNeighborhoodIterator<RealVectorImageType>

Type of the iterator that will be used to move through the image. Also the type which will be passed to the evaluate function

Definition at line 156 of file itkDisplacementFieldJacobianDeterminantFilter.h.

◆ ConstPointer

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::ConstPointer = SmartPointer<const Self>

◆ ImageBaseType

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::ImageBaseType = typename InputImageType::Superclass
protected

◆ InputImagePointer

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::InputImagePointer = typename InputImageType::Pointer

◆ InputImageType

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::InputImageType = TInputImage

Image type alias support

Definition at line 138 of file itkDisplacementFieldJacobianDeterminantFilter.h.

◆ InputPixelType

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::InputPixelType = typename TInputImage::PixelType

◆ OutputImagePointer

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::OutputImagePointer = typename OutputImageType::Pointer

◆ OutputImageType

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::OutputImageType = TOutputImage

◆ OutputPixelType

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::OutputPixelType = typename TOutputImage::PixelType

Extract some information from the image types. Dimensionality of the two images is assumed to be the same.

Definition at line 134 of file itkDisplacementFieldJacobianDeterminantFilter.h.

◆ Pointer

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::Pointer = SmartPointer<Self>

◆ RadiusType

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::RadiusType = typename ConstNeighborhoodIteratorType::RadiusType

◆ RealType

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::RealType = TRealType

Define the data type and the vector of data type used in calculations.

Definition at line 150 of file itkDisplacementFieldJacobianDeterminantFilter.h.

◆ RealVectorImageType

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::RealVectorImageType = Image<RealVectorType, TInputImage::ImageDimension>

◆ RealVectorType

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::RealVectorType = Vector<TRealType, InputPixelType::Dimension>

◆ Self

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::Self = DisplacementFieldJacobianDeterminantFilter

Standard class type aliases.

Definition at line 121 of file itkDisplacementFieldJacobianDeterminantFilter.h.

◆ Superclass

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::Superclass = ImageToImageFilter<TInputImage, TOutputImage>

◆ WeightsType

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
using itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::WeightsType = FixedArray<TRealType, ImageDimension>

Constructor & Destructor Documentation

◆ DisplacementFieldJacobianDeterminantFilter()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::DisplacementFieldJacobianDeterminantFilter ( )
protected

◆ ~DisplacementFieldJacobianDeterminantFilter()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::~DisplacementFieldJacobianDeterminantFilter ( )
overrideprotecteddefault

Member Function Documentation

◆ BeforeThreadedGenerateData()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
void itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::BeforeThreadedGenerateData ( )
overrideprotectedvirtual

Do any necessary casting/copying of the input data. Input pixel types whose value types are not real number types must be cast to real number types.

Reimplemented from itk::ImageSource< TOutputImage >.

◆ DynamicThreadedGenerateData()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
void itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::DynamicThreadedGenerateData ( const OutputImageRegionType outputRegionForThread)
overrideprotectedvirtual

DisplacementFieldJacobianDeterminantFilter can be implemented as a multithreaded filter (we're only using vnl_det(), which is trivially thread safe). Therefore, this implementation provides a DynamicThreadedGenerateData() routine which is called for each processing thread. The output image data is allocated automatically by the superclass prior to calling DynamicThreadedGenerateData(). DynamicThreadedGenerateData can only write to the portion of the output image specified by the parameter "outputRegionForThread"

See also
ImageToImageFilter::ThreadedGenerateData(), ImageToImageFilter::GenerateData()

Reimplemented from itk::ImageSource< TOutputImage >.

◆ EvaluateAtNeighborhood()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
virtual TRealType itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::EvaluateAtNeighborhood ( const ConstNeighborhoodIteratorType it) const
protectedvirtual

◆ GenerateInputRequestedRegion()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
void itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::GenerateInputRequestedRegion ( )
overridevirtual

DisplacementFieldJacobianDeterminantFilter needs a larger input requested region than the output requested region (larger by the kernel size to calculate derivatives). As such, DisplacementFieldJacobianDeterminantFilter needs to provide an implementation for GenerateInputRequestedRegion() in order to inform the pipeline execution model.

See also
ImageToImageFilter::GenerateInputRequestedRegion()

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

◆ GetDerivativeWeights()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
virtual const WeightsType& itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::GetDerivativeWeights ( ) const
virtual

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

◆ GetNameOfClass()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
const char* itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::GetNameOfClass ( ) const
overridevirtual

◆ GetNeighborhoodRadius()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
virtual const RadiusType& itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::GetNeighborhoodRadius ( ) const
protectedvirtual

Get/Set the neighborhood radius used for gradient computation

◆ GetRealValuedInputImage()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
virtual const ImageBaseType* itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::GetRealValuedInputImage ( ) const
protectedvirtual

Get access to the input image casted as real pixel values

◆ GetUseImageSpacing()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
virtual bool itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::GetUseImageSpacing ( ) const
virtual

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.

◆ New()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
static Pointer itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::New ( )
static

Method for creation through the object factory.

◆ PrintSelf()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
void itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::PrintSelf ( std::ostream &  os,
Indent  indent 
) const
overrideprotectedvirtual

Methods invoked by Print() to print information about the object including superclasses. Typically not called by the user (use Print() instead) but used in the hierarchical print process to combine the output of several classes.

Reimplemented from itk::ImageToImageFilter< TInputImage, TOutputImage >.

◆ SetDerivativeWeights()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
void itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::SetDerivativeWeights ( const WeightsType )

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

◆ SetNeighborhoodRadius()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
virtual void itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::SetNeighborhoodRadius ( RadiusType  _arg)
protectedvirtual

Get/Set the neighborhood radius used for gradient computation

◆ SetUseImageSpacing()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
void itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::SetUseImageSpacing ( bool  )

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()

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
virtual void itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::UseImageSpacingOn ( )
virtual

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.

Member Data Documentation

◆ ImageDimension

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
constexpr unsigned int itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::ImageDimension = TOutputImage::ImageDimension
staticconstexpr

The dimensionality of the input and output images.

Definition at line 144 of file itkDisplacementFieldJacobianDeterminantFilter.h.

◆ m_DerivativeWeights

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
WeightsType itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::m_DerivativeWeights {}
protected

The weights used to scale partial derivatives during processing

Definition at line 260 of file itkDisplacementFieldJacobianDeterminantFilter.h.

◆ m_HalfDerivativeWeights

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
WeightsType itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::m_HalfDerivativeWeights {}
protected

Pre-compute 0.5*m_DerivativeWeights since that is the only thing used in the computations.

Definition at line 264 of file itkDisplacementFieldJacobianDeterminantFilter.h.

◆ m_NeighborhoodRadius

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
RadiusType itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::m_NeighborhoodRadius {}
private

◆ m_RealValuedInputImage

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
ImageBaseType::ConstPointer itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::m_RealValuedInputImage {}
private

◆ m_RequestedNumberOfWorkUnits

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
ThreadIdType itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::m_RequestedNumberOfWorkUnits {}
private

◆ m_UseImageSpacing

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
bool itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::m_UseImageSpacing { true }
private

◆ VectorDimension

template<typename TInputImage , typename TRealType = float, typename TOutputImage = Image<TRealType, TInputImage::ImageDimension>>
constexpr unsigned int itk::DisplacementFieldJacobianDeterminantFilter< TInputImage, TRealType, TOutputImage >::VectorDimension = InputPixelType::Dimension
staticconstexpr

Length of the vector pixel type of the input image.

Definition at line 147 of file itkDisplacementFieldJacobianDeterminantFilter.h.


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