ITK  4.2.0
Insight Segmentation and Registration Toolkit
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType > Class Template Reference

#include <itkBSplineInterpolateImageFunction.h>

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

List of all members.

Public Types

typedef TCoefficientType CoefficientDataType
typedef
BSplineDecompositionImageFilter
< TImageType,
CoefficientImageType
CoefficientFilter
typedef CoefficientFilter::Pointer CoefficientFilterPointer
typedef Image
< CoefficientDataType,
itkGetStaticConstMacro(ImageDimension) > 
CoefficientImageType
typedef SmartPointer< const SelfConstPointer
typedef
Superclass::ContinuousIndexType 
ContinuousIndexType
typedef CovariantVector
< OutputType,
itkGetStaticConstMacro(ImageDimension) > 
CovariantVectorType
typedef Superclass::IndexType IndexType
typedef Superclass::InputImageType InputImageType
typedef
ImageLinearIteratorWithIndex
< TImageType > 
Iterator
typedef Superclass::OutputType OutputType
typedef SmartPointer< SelfPointer
typedef Superclass::PointType PointType
typedef
BSplineInterpolateImageFunction 
Self
typedef
InterpolateImageFunction
< TImageType, TCoordRep > 
Superclass
- Public Types inherited from itk::InterpolateImageFunction< TImageType, TCoordRep >
typedef Superclass::IndexValueType IndexValueType
typedef NumericTraits
< typename
TImageType::PixelType >
::RealType 
RealType
- Public Types inherited from itk::ImageFunction< TImageType, NumericTraits< TImageType::PixelType >::RealType, TCoordRep >
typedef TCoordRep CoordRepType
typedef
InputImageType::ConstPointer 
InputImageConstPointer
typedef InputImageType::PixelType InputPixelType
- Public Types inherited from itk::FunctionBase< Point< TCoordRep, TImageType::ImageDimension >, NumericTraits< TImageType::PixelType >::RealType >
typedef Point< TCoordRep,
TImageType::ImageDimension > 
InputType
- Public Types inherited from itk::Object
- Public Types inherited from itk::LightObject

Public Member Functions

virtual ::itk::LightObject::Pointer CreateAnother (void) const
virtual OutputType Evaluate (const PointType &point) const
virtual OutputType Evaluate (const PointType &point, ThreadIdType threadID) const
virtual OutputType EvaluateAtContinuousIndex (const ContinuousIndexType &index) const
virtual OutputType EvaluateAtContinuousIndex (const ContinuousIndexType &index, 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 &deriv, ThreadIdType threadID) const
virtual const char * GetNameOfClass () const
virtual ThreadIdType GetNumberOfThreads () const
virtual int GetSplineOrder () const
virtual void SetInputImage (const TImageType *inputData)
void SetNumberOfThreads (ThreadIdType numThreads)
void SetSplineOrder (unsigned int SplineOrder)
virtual void SetUseImageDirection (bool _arg)
virtual bool GetUseImageDirection () const
virtual void UseImageDirectionOn ()
virtual void UseImageDirectionOff ()
- Public Member Functions inherited from itk::InterpolateImageFunction< TImageType, TCoordRep >
virtual OutputType Evaluate (const PointType &point) const
virtual OutputType EvaluateAtContinuousIndex (const ContinuousIndexType &index) const =0
virtual OutputType EvaluateAtIndex (const IndexType &index) const
- Public Member Functions inherited from itk::ImageFunction< TImageType, NumericTraits< TImageType::PixelType >::RealType, TCoordRep >
void ConvertContinuousIndexToNearestIndex (const ContinuousIndexType &cindex, IndexType &index) const
void ConvertPointToContinuousIndex (const PointType &point, ContinuousIndexType &cindex) const
void ConvertPointToNearestIndex (const PointType &point, IndexType &index) const
virtual const ContinuousIndexTypeGetEndContinuousIndex ()
virtual const IndexTypeGetEndIndex ()
const InputImageTypeGetInputImage () const
virtual const ContinuousIndexTypeGetStartContinuousIndex ()
virtual const IndexTypeGetStartIndex ()
virtual bool IsInsideBuffer (const IndexType &index) const
virtual bool IsInsideBuffer (const ContinuousIndexType &index) const
virtual bool IsInsideBuffer (const PointType &point) const

Static Public Member Functions

static Pointer New ()

Static Public Attributes

static const unsigned int ImageDimension = Superclass::ImageDimension
- Static Public Attributes inherited from itk::InterpolateImageFunction< TImageType, TCoordRep >
- Static Public Attributes inherited from itk::ImageFunction< TImageType, NumericTraits< TImageType::PixelType >::RealType, TCoordRep >

Protected Member Functions

 BSplineInterpolateImageFunction ()
virtual OutputType EvaluateAtContinuousIndexInternal (const ContinuousIndexType &index, 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
 ~BSplineInterpolateImageFunction ()
- Protected Member Functions inherited from itk::InterpolateImageFunction< TImageType, TCoordRep >
 InterpolateImageFunction ()
 ~InterpolateImageFunction ()
- Protected Member Functions inherited from itk::ImageFunction< TImageType, NumericTraits< TImageType::PixelType >::RealType, TCoordRep >
 ImageFunction ()
 ~ImageFunction ()
- Protected Member Functions inherited from itk::FunctionBase< Point< TCoordRep, TImageType::ImageDimension >, NumericTraits< TImageType::PixelType >::RealType >
 FunctionBase ()
 ~FunctionBase ()
- Protected Member Functions inherited from itk::Object
 Object ()
bool PrintObservers (std::ostream &os, Indent indent) const
virtual void SetTimeStamp (const TimeStamp &time)
virtual ~Object ()
- 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

Private Member Functions

void ApplyMirrorBoundaryConditions (vnl_matrix< long > &evaluateIndex, unsigned int splineOrder) const
 BSplineInterpolateImageFunction (const Self &)
void DetermineRegionOfSupport (vnl_matrix< long > &evaluateIndex, const ContinuousIndexType &x, unsigned int splineOrder) const
void GeneratePointsToIndex ()
void operator= (const Self &)
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_NumberOfThreads
std::vector< IndexTypem_PointsToIndex
vnl_matrix< long > * m_ThreadedEvaluateIndex
vnl_matrix< double > * m_ThreadedWeights
vnl_matrix< double > * m_ThreadedWeightsDerivative
bool m_UseImageDirection

Detailed Description

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
class itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, 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
Wiki Examples:

Definition at line 84 of file itkBSplineInterpolateImageFunction.h.


Member Typedef Documentation

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef TCoefficientType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::CoefficientDataType

Internal Coefficient typedef support

Definition at line 122 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef BSplineDecompositionImageFilter< TImageType, CoefficientImageType > itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::CoefficientFilter

Define filter for calculating the BSpline coefficients

Definition at line 128 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef CoefficientFilter::Pointer itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::CoefficientFilterPointer

Definition at line 129 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef Image< CoefficientDataType, itkGetStaticConstMacro(ImageDimension) > itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::CoefficientImageType

Definition at line 125 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef SmartPointer< const Self > itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::ConstPointer
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef Superclass::ContinuousIndexType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::ContinuousIndexType
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef CovariantVector< OutputType, itkGetStaticConstMacro(ImageDimension) > itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::CovariantVectorType

Derivative typedef support

Definition at line 134 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef Superclass::IndexType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::IndexType

Index typedef support.

Reimplemented from itk::InterpolateImageFunction< TImageType, TCoordRep >.

Definition at line 110 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef Superclass::InputImageType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::InputImageType

InputImageType typedef support.

Reimplemented from itk::InterpolateImageFunction< TImageType, TCoordRep >.

Definition at line 104 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef ImageLinearIteratorWithIndex< TImageType > itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::Iterator

Iterator typedef support

Definition at line 119 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef Superclass::OutputType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::OutputType

OutputType typedef support.

Reimplemented from itk::InterpolateImageFunction< TImageType, TCoordRep >.

Definition at line 98 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef SmartPointer< Self > itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::Pointer
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef Superclass::PointType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::PointType

PointType typedef support

Reimplemented from itk::InterpolateImageFunction< TImageType, TCoordRep >.

Definition at line 116 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef BSplineInterpolateImageFunction itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::Self
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
typedef InterpolateImageFunction< TImageType, TCoordRep > itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::Superclass

Constructor & Destructor Documentation

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::BSplineInterpolateImageFunction ( )
protected
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::~BSplineInterpolateImageFunction ( )
protected
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::BSplineInterpolateImageFunction ( const Self )
private

Member Function Documentation

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

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

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual::itk::LightObject::Pointer itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::CreateAnother ( void  ) const
virtual

Create an object from an instance, potentially deferring to a factory. This method allows you to create an instance of an object that is exactly the same type as the referring object. This is useful in cases where an object has been cast back to a base class.

Reimplemented from itk::Object.

Reimplemented in itk::BSplineResampleImageFunction< TImageType, TCoordRep >.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, 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

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual OutputType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::Evaluate ( const PointType point) const
inlinevirtual

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.

Definition at line 144 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual OutputType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::Evaluate ( const PointType point,
ThreadIdType  threadID 
) const
inlinevirtual

Definition at line 154 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual OutputType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::EvaluateAtContinuousIndex ( const ContinuousIndexType index) const
inlinevirtual

Definition at line 164 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual OutputType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::EvaluateAtContinuousIndex ( const ContinuousIndexType index,
ThreadIdType  threadID 
) const
virtual
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual OutputType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::EvaluateAtContinuousIndexInternal ( const ContinuousIndexType index,
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.

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

Definition at line 183 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
CovariantVectorType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::EvaluateDerivative ( const PointType point,
ThreadIdType  threadID 
) const
inline

Definition at line 193 of file itkBSplineInterpolateImageFunction.h.

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

Definition at line 203 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
CovariantVectorType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::EvaluateDerivativeAtContinuousIndex ( const ContinuousIndexType x,
ThreadIdType  threadID 
) const
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual CovariantVectorType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::EvaluateDerivativeAtContinuousIndexInternal ( const ContinuousIndexType x,
vnl_matrix< long > &  evaluateIndex,
vnl_matrix< double > &  weights,
vnl_matrix< double > &  weightsDerivative 
) const
protectedvirtual
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::EvaluateValueAndDerivative ( const PointType point,
OutputType value,
CovariantVectorType deriv 
) const
inline

Definition at line 227 of file itkBSplineInterpolateImageFunction.h.

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

Definition at line 242 of file itkBSplineInterpolateImageFunction.h.

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

Definition at line 257 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::EvaluateValueAndDerivativeAtContinuousIndex ( const ContinuousIndexType x,
OutputType value,
CovariantVectorType deriv,
ThreadIdType  threadID 
) const
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::EvaluateValueAndDerivativeAtContinuousIndexInternal ( const ContinuousIndexType x,
OutputType value,
CovariantVectorType derivativeValue,
vnl_matrix< long > &  evaluateIndex,
vnl_matrix< double > &  weights,
vnl_matrix< double > &  weightsDerivative 
) const
protectedvirtual
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::GeneratePointsToIndex ( )
private

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

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual const char* itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::GetNameOfClass ( ) const
virtual

Run-time type information (and related methods).

Reimplemented from itk::InterpolateImageFunction< TImageType, TCoordRep >.

Reimplemented in itk::BSplineResampleImageFunction< TImageType, TCoordRep >.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual ThreadIdType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::GetNumberOfThreads ( ) const
virtual
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual int itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::GetSplineOrder ( ) const
virtual
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual bool itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, 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 coodinate 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.

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

New macro for creation of through a Smart Pointer

Reimplemented from itk::Object.

Reimplemented in itk::BSplineResampleImageFunction< TImageType, TCoordRep >.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::operator= ( const Self )
private

Mutex lock to protect modification to the reference count

Reimplemented from itk::InterpolateImageFunction< TImageType, TCoordRep >.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::PrintSelf ( std::ostream &  os,
Indent  indent 
) const
protectedvirtual

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::InterpolateImageFunction< TImageType, TCoordRep >.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, 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

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::SetInputImage ( const TImageType *  inputData)
virtual
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, 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

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::SetNumberOfThreads ( ThreadIdType  numThreads)
template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::SetSplineOrder ( unsigned int  SplineOrder)

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

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, 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 coodinate 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.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::UseImageDirectionOff ( )
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 coodinate 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.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
virtual void itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, 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 coodinate 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

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
const unsigned int itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::ImageDimension = Superclass::ImageDimension
static

Dimension underlying input image.

Reimplemented from itk::InterpolateImageFunction< TImageType, TCoordRep >.

Definition at line 107 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
Iterator itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_CIterator
private

Definition at line 397 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
CoefficientFilterPointer itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_CoefficientFilter
private

Definition at line 409 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
CoefficientImageType::ConstPointer itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_Coefficients
protected

Definition at line 366 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
TImageType::SizeType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_DataLength
protected

Definition at line 361 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
unsigned long itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_MaxNumberInterpolationPoints
private

Definition at line 400 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
ThreadIdType itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_NumberOfThreads
private

Definition at line 415 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
std::vector< IndexType > itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_PointsToIndex
private

Definition at line 404 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
std::vector< CoefficientDataType > itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_Scratch
protected

Definition at line 359 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
unsigned int itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_SplineOrder
protected

Definition at line 363 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
vnl_matrix< long >* itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_ThreadedEvaluateIndex
private

Definition at line 416 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
vnl_matrix< double >* itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_ThreadedWeights
private

Definition at line 417 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
vnl_matrix< double >* itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_ThreadedWeightsDerivative
private

Definition at line 418 of file itkBSplineInterpolateImageFunction.h.

template<class TImageType, class TCoordRep = double, class TCoefficientType = double>
bool itk::BSplineInterpolateImageFunction< TImageType, TCoordRep, TCoefficientType >::m_UseImageDirection
private

Definition at line 413 of file itkBSplineInterpolateImageFunction.h.


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