ITK  6.0.0
Insight Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType > Class Template Reference

#include <itkBSplineInterpolateImageFunction.h>

Detailed Description

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
class itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >

Evaluates the B-Spline interpolation of an image. Spline order may be from 0 to 5.

This class defines N-Dimension B-Spline transformation. It is based on:

[1] M. Unser,
    "Splines: A Perfect Fit for Signal and Image Processing,"
    IEEE Signal Processing Magazine, vol. 16, no. 6, pp. 22-38,
    November 1999.
[2] M. Unser, A. Aldroubi and M. Eden,
    "B-Spline Signal Processing: Part I--Theory,"
    IEEE Transactions on Signal Processing, vol. 41, no. 2, pp. 821-832,
    February 1993.
[3] M. Unser, A. Aldroubi and M. Eden,
    "B-Spline Signal Processing: Part II--Efficient Design and Applications,"
     IEEE Transactions on Signal Processing, vol. 41, no. 2, pp. 834-848,
     February 1993.

And code obtained from bigwww.epfl.ch by Philippe Thevenaz

The B spline coefficients are calculated through the BSplineDecompositionImageFilter

Limitations: Spline order must be between 0 and 5. Spline order must be set before setting the image. Uses mirror boundary conditions. Requires the same order of Spline for each dimension. Spline is determined in all dimensions, cannot selectively pick dimension for calculating spline.

See also
BSplineDecompositionImageFilter
Examples
Examples/Filtering/ResampleImageFilter7.cxx, and SphinxExamples/src/Filtering/ImageGrid/UpsampleAnImage/Code.cxx.

Definition at line 83 of file itkBSplineInterpolateImageFunction.h.

+ Inheritance diagram for itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >:
+ Collaboration diagram for itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >:

Public Types

using CoefficientDataType = TCoefficientType
 
using CoefficientFilter = BSplineDecompositionImageFilter< TImageType, CoefficientImageType >
 
using CoefficientFilterPointer = typename CoefficientFilter::Pointer
 
using CoefficientImageType = Image< CoefficientDataType, Self::ImageDimension >
 
using ConstPointer = SmartPointer< const Self >
 
using CovariantVectorType = CovariantVector< OutputType, Self::ImageDimension >
 
using Iterator = ImageLinearIteratorWithIndex< TImageType >
 
using Pointer = SmartPointer< Self >
 
using Self = BSplineInterpolateImageFunction
 
using Superclass = InterpolateImageFunction< TImageType, TCoordinate >
 
- Public Types inherited from itk::InterpolateImageFunction< TImageType, TCoordinate >
using ConstPointer = SmartPointer< const Self >
 
using Pointer = SmartPointer< Self >
 
using RealType = typename NumericTraits< typename TImageType ::PixelType >::RealType
 
using Self = InterpolateImageFunction
 
using SizeType = typename InputImageType::SizeType
 
using Superclass = ImageFunction< TImageType, typename NumericTraits< typename TImageType ::PixelType >::RealType, TCoordinate >
 
- Public Types inherited from itk::ImageFunction< TImageType, NumericTraits< TImageType ::PixelType >::RealType, TCoordinate >
using ConstPointer = SmartPointer< const Self >
 
using ContinuousIndexType = ContinuousIndex< TCoordinate, Self::ImageDimension >
 
using CoordinateType = TCoordinate
 
using IndexType = typename InputImageType::IndexType
 
using IndexValueType = typename InputImageType::IndexValueType
 
using InputImageConstPointer = typename InputImageType::ConstPointer
 
using InputImageType = TImageType
 
using InputPixelType = typename InputImageType::PixelType
 
using OutputType = NumericTraits< TImageType ::PixelType >::RealType
 
using Pointer = SmartPointer< Self >
 
using PointType = Point< TCoordinate, Self::ImageDimension >
 
using Self = ImageFunction
 
using Superclass = FunctionBase< Point< TCoordinate, Self::ImageDimension >, NumericTraits< TImageType ::PixelType >::RealType >
 
- Public Types inherited from itk::FunctionBase< Point< TCoordinate, TImageType ::ImageDimension >, NumericTraits< TImageType ::PixelType >::RealType >
using ConstPointer = SmartPointer< const Self >
 
using InputType = Point< TCoordinate, TImageType ::ImageDimension >
 
using OutputType = NumericTraits< TImageType ::PixelType >::RealType
 
using Pointer = SmartPointer< Self >
 
using Self = FunctionBase
 
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

virtual OutputType Evaluate (const PointType &point, ThreadIdType threadId) const
 
OutputType EvaluateAtContinuousIndex (const ContinuousIndexType &index) const override
 
virtual OutputType EvaluateAtContinuousIndex (const ContinuousIndexType &x, ThreadIdType threadId) const
 
CovariantVectorType EvaluateDerivative (const PointType &point) const
 
CovariantVectorType EvaluateDerivative (const PointType &point, ThreadIdType threadId) const
 
CovariantVectorType EvaluateDerivativeAtContinuousIndex (const ContinuousIndexType &x) const
 
CovariantVectorType EvaluateDerivativeAtContinuousIndex (const ContinuousIndexType &x, ThreadIdType threadId) const
 
void EvaluateValueAndDerivative (const PointType &point, OutputType &value, CovariantVectorType &deriv) const
 
void EvaluateValueAndDerivative (const PointType &point, OutputType &value, CovariantVectorType &deriv, ThreadIdType threadId) const
 
void EvaluateValueAndDerivativeAtContinuousIndex (const ContinuousIndexType &x, OutputType &value, CovariantVectorType &deriv) const
 
void EvaluateValueAndDerivativeAtContinuousIndex (const ContinuousIndexType &x, OutputType &value, CovariantVectorType &derivativeValue, ThreadIdType threadId) const
 
const char * GetNameOfClass () const override
 
virtual ThreadIdType GetNumberOfWorkUnits () const
 
SizeType GetRadius () const override
 
virtual unsigned int GetSplineOrder () const
 
void SetInputImage (const TImageType *inputData) override
 
void SetNumberOfWorkUnits (ThreadIdType numWorkUnits)
 
void SetSplineOrder (unsigned int SplineOrder)
 
OutputType Evaluate (const PointType &point) const override
 
virtual void SetUseImageDirection (bool _arg)
 
virtual bool GetUseImageDirection () const
 
virtual void UseImageDirectionOn ()
 
- Public Member Functions inherited from itk::InterpolateImageFunction< TImageType, TCoordinate >
OutputType EvaluateAtIndex (const IndexType &index) const override
 
const char * GetNameOfClass () const override
 
OutputType Evaluate (const PointType &point) const override
 
- Public Member Functions inherited from itk::ImageFunction< TImageType, NumericTraits< TImageType ::PixelType >::RealType, TCoordinate >
void ConvertContinuousIndexToNearestIndex (const ContinuousIndexType &cindex, IndexType &index) const
 
void ConvertPointToContinuousIndex (const PointType &point, ContinuousIndexType &cindex) const
 
virtual const ContinuousIndexTypeGetEndContinuousIndex () const
 
virtual const IndexTypeGetEndIndex () const
 
const InputImageTypeGetInputImage () const
 
virtual const ContinuousIndexTypeGetStartContinuousIndex () const
 
virtual const IndexTypeGetStartIndex () const
 
virtual bool IsInsideBuffer (const IndexType &index) const
 
virtual bool IsInsideBuffer (const ContinuousIndexType &index) const
 
virtual bool IsInsideBuffer (const PointType &point) const
 
void ConvertPointToNearestIndex (const PointType &point, IndexType &index) 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::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 = Superclass::ImageDimension
 
- Static Public Attributes inherited from itk::InterpolateImageFunction< TImageType, TCoordinate >
static constexpr unsigned int ImageDimension
 
- Static Public Attributes inherited from itk::ImageFunction< TImageType, NumericTraits< TImageType ::PixelType >::RealType, TCoordinate >
static constexpr unsigned int ImageDimension
 

Protected Member Functions

 BSplineInterpolateImageFunction ()
 
virtual OutputType EvaluateAtContinuousIndexInternal (const ContinuousIndexType &x, vnl_matrix< long > &evaluateIndex, vnl_matrix< double > &weights) const
 
virtual CovariantVectorType EvaluateDerivativeAtContinuousIndexInternal (const ContinuousIndexType &x, vnl_matrix< long > &evaluateIndex, vnl_matrix< double > &weights, vnl_matrix< double > &weightsDerivative) const
 
virtual void EvaluateValueAndDerivativeAtContinuousIndexInternal (const ContinuousIndexType &x, OutputType &value, CovariantVectorType &derivativeValue, vnl_matrix< long > &evaluateIndex, vnl_matrix< double > &weights, vnl_matrix< double > &weightsDerivative) const
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
 ~BSplineInterpolateImageFunction () override=default
 
- Protected Member Functions inherited from itk::InterpolateImageFunction< TImageType, TCoordinate >
 InterpolateImageFunction ()=default
 
void PrintSelf (std::ostream &os, Indent indent) const override
 
 ~InterpolateImageFunction () override=default
 
- Protected Member Functions inherited from itk::ImageFunction< TImageType, NumericTraits< TImageType ::PixelType >::RealType, TCoordinate >
 ImageFunction ()
 
 ~ImageFunction () override=default
 
- Protected Member Functions inherited from itk::FunctionBase< Point< TCoordinate, TImageType ::ImageDimension >, NumericTraits< TImageType ::PixelType >::RealType >
 FunctionBase ()=default
 
 ~FunctionBase () override=default
 
- 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

CoefficientImageType::ConstPointer m_Coefficients {}
 
TImageType::SizeType m_DataLength {}
 
std::vector< CoefficientDataTypem_Scratch {}
 
unsigned int m_SplineOrder {}
 
- Protected Attributes inherited from itk::ImageFunction< TImageType, NumericTraits< TImageType ::PixelType >::RealType, TCoordinate >
ContinuousIndexType m_EndContinuousIndex
 
IndexType m_EndIndex
 
InputImageConstPointer m_Image
 
ContinuousIndexType m_StartContinuousIndex
 
IndexType m_StartIndex
 
- Protected Attributes inherited from itk::LightObject
std::atomic< int > m_ReferenceCount {}
 

Private Member Functions

void ApplyMirrorBoundaryConditions (vnl_matrix< long > &evaluateIndex, unsigned int splineOrder) const
 
void DetermineRegionOfSupport (vnl_matrix< long > &evaluateIndex, const ContinuousIndexType &x, unsigned int splineOrder) const
 
void GeneratePointsToIndex ()
 
void SetDerivativeWeights (const ContinuousIndexType &x, const vnl_matrix< long > &EvaluateIndex, vnl_matrix< double > &weights, unsigned int splineOrder) const
 
void SetInterpolationWeights (const ContinuousIndexType &x, const vnl_matrix< long > &EvaluateIndex, vnl_matrix< double > &weights, unsigned int splineOrder) const
 

Private Attributes

Iterator m_CIterator {}
 
CoefficientFilterPointer m_CoefficientFilter {}
 
unsigned long m_MaxNumberInterpolationPoints {}
 
ThreadIdType m_NumberOfWorkUnits {}
 
std::vector< IndexTypem_PointsToIndex {}
 
std::unique_ptr< vnl_matrix< long >[]> m_ThreadedEvaluateIndex
 
std::unique_ptr< vnl_matrix< double >[]> m_ThreadedWeights
 
std::unique_ptr< vnl_matrix< double >[]> m_ThreadedWeightsDerivative
 
bool m_UseImageDirection { true }
 

Member Typedef Documentation

◆ CoefficientDataType

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
using itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::CoefficientDataType = TCoefficientType

Internal Coefficient type alias support

Definition at line 125 of file itkBSplineInterpolateImageFunction.h.

◆ CoefficientFilter

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
using itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::CoefficientFilter = BSplineDecompositionImageFilter<TImageType, CoefficientImageType>

Define filter for calculating the BSpline coefficients

Definition at line 129 of file itkBSplineInterpolateImageFunction.h.

◆ CoefficientFilterPointer

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
using itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::CoefficientFilterPointer = typename CoefficientFilter::Pointer

Definition at line 130 of file itkBSplineInterpolateImageFunction.h.

◆ CoefficientImageType

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
using itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::CoefficientImageType = Image<CoefficientDataType, Self::ImageDimension>

Definition at line 126 of file itkBSplineInterpolateImageFunction.h.

◆ ConstPointer

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
using itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::ConstPointer = SmartPointer<const Self>

Definition at line 92 of file itkBSplineInterpolateImageFunction.h.

◆ CovariantVectorType

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
using itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::CovariantVectorType = CovariantVector<OutputType, Self::ImageDimension>

Derivative type alias support

Definition at line 133 of file itkBSplineInterpolateImageFunction.h.

◆ Iterator

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
using itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::Iterator = ImageLinearIteratorWithIndex<TImageType>

Iterator type alias support

Definition at line 122 of file itkBSplineInterpolateImageFunction.h.

◆ Pointer

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
using itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::Pointer = SmartPointer<Self>

Definition at line 91 of file itkBSplineInterpolateImageFunction.h.

◆ Self

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
using itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::Self = BSplineInterpolateImageFunction

Standard class type aliases.

Definition at line 89 of file itkBSplineInterpolateImageFunction.h.

◆ Superclass

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
using itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::Superclass = InterpolateImageFunction<TImageType, TCoordinate>

Definition at line 90 of file itkBSplineInterpolateImageFunction.h.

Constructor & Destructor Documentation

◆ BSplineInterpolateImageFunction()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::BSplineInterpolateImageFunction ( )
protected

◆ ~BSplineInterpolateImageFunction()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::~BSplineInterpolateImageFunction ( )
overrideprotecteddefault

Member Function Documentation

◆ ApplyMirrorBoundaryConditions()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::ApplyMirrorBoundaryConditions ( vnl_matrix< long > &  evaluateIndex,
unsigned int  splineOrder 
) const
private

Set the indices in evaluateIndex at the boundaries based on mirror boundary conditions.

◆ DetermineRegionOfSupport()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::DetermineRegionOfSupport ( vnl_matrix< long > &  evaluateIndex,
const ContinuousIndexType x,
unsigned int  splineOrder 
) const
private

Determines the indices to use give the splines region of support

◆ Evaluate() [1/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
OutputType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::Evaluate ( const PointType point) const
inlineoverridevirtual

Evaluate the function at a ContinuousIndex position.

Returns the B-Spline interpolated image intensity at a
specified point position. No bounds checking is done.
The point is assume to lie within the image buffer.

ImageFunction::IsInsideBuffer() can be used to check bounds before
calling the method. 

Implements itk::ImageFunction< TImageType, NumericTraits< TImageType ::PixelType >::RealType, TCoordinate >.

Definition at line 144 of file itkBSplineInterpolateImageFunction.h.

◆ Evaluate() [2/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
virtual OutputType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::Evaluate ( const PointType point,
ThreadIdType  threadId 
) const
inlinevirtual

Definition at line 154 of file itkBSplineInterpolateImageFunction.h.

◆ EvaluateAtContinuousIndex() [1/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
OutputType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateAtContinuousIndex ( const ContinuousIndexType index) const
inlineoverridevirtual

Interpolate the image at a continuous index position

Returns the interpolated image intensity at a specified index position. No bounds checking is done. The point is assume to lie within the image buffer.

Subclasses must override this method.

ImageFunction::IsInsideBuffer() can be used to check bounds before calling the method.

Implements itk::InterpolateImageFunction< TImageType, TCoordinate >.

Definition at line 162 of file itkBSplineInterpolateImageFunction.h.

◆ EvaluateAtContinuousIndex() [2/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
virtual OutputType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateAtContinuousIndex ( const ContinuousIndexType x,
ThreadIdType  threadId 
) const
inlinevirtual

Definition at line 175 of file itkBSplineInterpolateImageFunction.h.

◆ EvaluateAtContinuousIndexInternal()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
virtual OutputType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateAtContinuousIndexInternal ( const ContinuousIndexType x,
vnl_matrix< long > &  evaluateIndex,
vnl_matrix< double > &  weights 
) const
protectedvirtual

The following methods take working space (evaluateIndex, weights, weightsDerivative) that is managed by the caller. If threadId is known, the working variables are looked up in the thread indexed arrays. If threadId is not known, working variables are made on the stack and passed to these methods. The stack allocation should be ok since these methods do not store the working variables, i.e. they are not expected to be available beyond the scope of the function call.

This was done to allow for two types of re-entrancy. The first is when a threaded filter, e.g. InterpolateImagePointsFilter calls EvaluateAtContinuousIndex from multiple threads without passing a threadId. So, EvaluateAtContinuousIndex must be thread safe. This is handled with the stack-based allocation of the working space.

The second form of re-entrancy involves methods that call EvaluateAtContinuousIndex from multiple threads, but pass a threadId. In this case, we can gain a little efficiency (hopefully) by looking up pre-allocated working space in arrays that are indexed by thread. The efficiency gain is likely dependent on the size of the working variables, which are in-turn dependent on the dimensionality of the image and the order of the spline.

◆ EvaluateDerivative() [1/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
CovariantVectorType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateDerivative ( const PointType point) const
inline

Definition at line 182 of file itkBSplineInterpolateImageFunction.h.

◆ EvaluateDerivative() [2/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
CovariantVectorType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateDerivative ( const PointType point,
ThreadIdType  threadId 
) const
inline

Definition at line 192 of file itkBSplineInterpolateImageFunction.h.

◆ EvaluateDerivativeAtContinuousIndex() [1/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
CovariantVectorType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateDerivativeAtContinuousIndex ( const ContinuousIndexType x) const
inline

Definition at line 200 of file itkBSplineInterpolateImageFunction.h.

◆ EvaluateDerivativeAtContinuousIndex() [2/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
CovariantVectorType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateDerivativeAtContinuousIndex ( const ContinuousIndexType x,
ThreadIdType  threadId 
) const
inline

Definition at line 217 of file itkBSplineInterpolateImageFunction.h.

◆ EvaluateDerivativeAtContinuousIndexInternal()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
virtual CovariantVectorType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateDerivativeAtContinuousIndexInternal ( const ContinuousIndexType x,
vnl_matrix< long > &  evaluateIndex,
vnl_matrix< double > &  weights,
vnl_matrix< double > &  weightsDerivative 
) const
protectedvirtual

◆ EvaluateValueAndDerivative() [1/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateValueAndDerivative ( const PointType point,
OutputType value,
CovariantVectorType deriv 
) const
inline

Definition at line 224 of file itkBSplineInterpolateImageFunction.h.

◆ EvaluateValueAndDerivative() [2/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateValueAndDerivative ( const PointType point,
OutputType value,
CovariantVectorType deriv,
ThreadIdType  threadId 
) const
inline

Definition at line 234 of file itkBSplineInterpolateImageFunction.h.

◆ EvaluateValueAndDerivativeAtContinuousIndex() [1/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateValueAndDerivativeAtContinuousIndex ( const ContinuousIndexType x,
OutputType value,
CovariantVectorType deriv 
) const
inline

Definition at line 245 of file itkBSplineInterpolateImageFunction.h.

◆ EvaluateValueAndDerivativeAtContinuousIndex() [2/2]

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateValueAndDerivativeAtContinuousIndex ( const ContinuousIndexType x,
OutputType value,
CovariantVectorType derivativeValue,
ThreadIdType  threadId 
) const
inline

Definition at line 265 of file itkBSplineInterpolateImageFunction.h.

◆ EvaluateValueAndDerivativeAtContinuousIndexInternal()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
virtual void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::EvaluateValueAndDerivativeAtContinuousIndexInternal ( const ContinuousIndexType x,
OutputType value,
CovariantVectorType derivativeValue,
vnl_matrix< long > &  evaluateIndex,
vnl_matrix< double > &  weights,
vnl_matrix< double > &  weightsDerivative 
) const
protectedvirtual

◆ GeneratePointsToIndex()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::GeneratePointsToIndex ( )
private

Precomputation for converting the 1D index of the interpolation neighborhood to an N-dimensional index.

◆ GetNameOfClass()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
const char* itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::GetNameOfClass ( ) const
overridevirtual
See also
LightObject::GetNameOfClass()

Reimplemented from itk::Object.

◆ GetNumberOfWorkUnits()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
virtual ThreadIdType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::GetNumberOfWorkUnits ( ) const
virtual

◆ GetRadius()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
SizeType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::GetRadius ( ) const
inlineoverridevirtual

Get the radius required for interpolation.

This defines the number of surrounding pixels required to interpolate at a given point.

Implements itk::InterpolateImageFunction< TImageType, TCoordinate >.

Definition at line 310 of file itkBSplineInterpolateImageFunction.h.

◆ GetSplineOrder()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
virtual unsigned int itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::GetSplineOrder ( ) const
virtual

◆ GetUseImageDirection()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
virtual bool itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::GetUseImageDirection ( ) const
virtual

The UseImageDirection flag determines whether image derivatives are computed with respect to the image grid or with respect to the physical space. When this flag is ON the derivatives are computed with respect to the coordinate system of physical space. The difference is whether we take into account the image Direction or not. The flag ON will take into account the image direction and will result in an extra matrix multiplication compared to the amount of computation performed when the flag is OFF. The default value of this flag is On.

◆ New()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
static Pointer itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::New ( )
static

New macro for creation of through a Smart Pointer

◆ PrintSelf()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::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::Object.

◆ SetDerivativeWeights()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::SetDerivativeWeights ( const ContinuousIndexType x,
const vnl_matrix< long > &  EvaluateIndex,
vnl_matrix< double > &  weights,
unsigned int  splineOrder 
) const
private

Determines the weights for the derivative portion of the value x

◆ SetInputImage()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::SetInputImage ( const TImageType *  inputData)
overridevirtual

Set the input image. This must be set by the user.

Reimplemented from itk::ImageFunction< TImageType, NumericTraits< TImageType ::PixelType >::RealType, TCoordinate >.

◆ SetInterpolationWeights()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::SetInterpolationWeights ( const ContinuousIndexType x,
const vnl_matrix< long > &  EvaluateIndex,
vnl_matrix< double > &  weights,
unsigned int  splineOrder 
) const
private

Determines the weights for interpolation of the value x

◆ SetNumberOfWorkUnits()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::SetNumberOfWorkUnits ( ThreadIdType  numWorkUnits)

◆ SetSplineOrder()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::SetSplineOrder ( unsigned int  SplineOrder)

Get/Sets the Spline Order, supports 0th - 5th order splines. The default is a 3rd order spline.

◆ SetUseImageDirection()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
virtual void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::SetUseImageDirection ( bool  _arg)
virtual

The UseImageDirection flag determines whether image derivatives are computed with respect to the image grid or with respect to the physical space. When this flag is ON the derivatives are computed with respect to the coordinate system of physical space. The difference is whether we take into account the image Direction or not. The flag ON will take into account the image direction and will result in an extra matrix multiplication compared to the amount of computation performed when the flag is OFF. The default value of this flag is On.

◆ UseImageDirectionOn()

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
virtual void itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::UseImageDirectionOn ( )
virtual

The UseImageDirection flag determines whether image derivatives are computed with respect to the image grid or with respect to the physical space. When this flag is ON the derivatives are computed with respect to the coordinate system of physical space. The difference is whether we take into account the image Direction or not. The flag ON will take into account the image direction and will result in an extra matrix multiplication compared to the amount of computation performed when the flag is OFF. The default value of this flag is On.

Member Data Documentation

◆ ImageDimension

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
constexpr unsigned int itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::ImageDimension = Superclass::ImageDimension
staticconstexpr

Dimension underlying input image.

Definition at line 107 of file itkBSplineInterpolateImageFunction.h.

◆ m_CIterator

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
Iterator itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_CIterator {}
private

Definition at line 400 of file itkBSplineInterpolateImageFunction.h.

◆ m_CoefficientFilter

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
CoefficientFilterPointer itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_CoefficientFilter {}
private

Definition at line 412 of file itkBSplineInterpolateImageFunction.h.

◆ m_Coefficients

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
CoefficientImageType::ConstPointer itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_Coefficients {}
protected

Definition at line 367 of file itkBSplineInterpolateImageFunction.h.

◆ m_DataLength

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
TImageType::SizeType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_DataLength {}
protected

Definition at line 362 of file itkBSplineInterpolateImageFunction.h.

◆ m_MaxNumberInterpolationPoints

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
unsigned long itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_MaxNumberInterpolationPoints {}
private

Definition at line 403 of file itkBSplineInterpolateImageFunction.h.

◆ m_NumberOfWorkUnits

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
ThreadIdType itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_NumberOfWorkUnits {}
private

Definition at line 418 of file itkBSplineInterpolateImageFunction.h.

◆ m_PointsToIndex

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
std::vector<IndexType> itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_PointsToIndex {}
private

Definition at line 407 of file itkBSplineInterpolateImageFunction.h.

◆ m_Scratch

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
std::vector<CoefficientDataType> itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_Scratch {}
protected

Definition at line 360 of file itkBSplineInterpolateImageFunction.h.

◆ m_SplineOrder

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
unsigned int itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_SplineOrder {}
protected

Definition at line 364 of file itkBSplineInterpolateImageFunction.h.

◆ m_ThreadedEvaluateIndex

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
std::unique_ptr<vnl_matrix<long>[]> itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_ThreadedEvaluateIndex
private

Definition at line 419 of file itkBSplineInterpolateImageFunction.h.

◆ m_ThreadedWeights

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
std::unique_ptr<vnl_matrix<double>[]> itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_ThreadedWeights
private

Definition at line 420 of file itkBSplineInterpolateImageFunction.h.

◆ m_ThreadedWeightsDerivative

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
std::unique_ptr<vnl_matrix<double>[]> itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_ThreadedWeightsDerivative
private

Definition at line 421 of file itkBSplineInterpolateImageFunction.h.

◆ m_UseImageDirection

template<typename TImageType, typename TCoordinate = double, typename TCoefficientType = double>
bool itk::BSplineInterpolateImageFunction< TImageType, TCoordinate, TCoefficientType >::m_UseImageDirection { true }
private

Definition at line 416 of file itkBSplineInterpolateImageFunction.h.


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