aimsalgo  5.0.5
Neuroimaging image processing
aims::SplineResampler< T > Class Template Referenceabstract

B-Spline-based resampling. More...

#include <aims/resampling/splineresampler.h>

Inheritance diagram for aims::SplineResampler< T >:
Collaboration diagram for aims::SplineResampler< T >:

Public Types

typedef carto::DataTypeTraits< T >::ChannelType ChannelType
 
- Public Types inherited from carto::RCObject
typedef int RefCounterType
 

Public Member Functions

 SplineResampler ()
 
 ~SplineResampler ()
 
virtual int getOrder () const =0
 Spline order (1 to 7) More...
 
carto::VolumeRef< double > getSplineCoef (const carto::Volume< T > &inVolume, int t=0, bool verbose=false)
 Computes spline coefficients corresponding to an input volume. More...
 
void reset ()
 Clear the cache. More...
 
void resample_inv_to_vox (const carto::Volume< T > &input_data, const soma::Transformation3d &inverse_transform_to_vox, const T &background, carto::Volume< T > &output_data, bool verbose=false) const CARTO_OVERRIDE
 Resample a volume into an existing output volume. More...
 
void resample_channel_inv_to_vox (const carto::Volume< ChannelType > &inVolume, const soma::Transformation3d &inverse_transform_to_vox, const ChannelType &outBackground, carto::Volume< ChannelType > &outVolume, bool verbose) const
 
- Public Member Functions inherited from aims::Resampler< T >
 Resampler ()
 
virtual ~Resampler ()
 
void doit (const aims::AffineTransformation3d &transform, carto::Volume< T > &output_data) const
 Resample the input volume set with setRef() into an existing volume. More...
 
carto::VolumeRef< T > doit (const aims::AffineTransformation3d &transform, int dimX, int dimY, int dimZ, const Point3df &voxel_size) const
 Resample the input volume set with setRef() in a newly allocated volume. More...
 
virtual void resample (const carto::Volume< T > &input_data, const aims::AffineTransformation3d &transform, const T &background, carto::Volume< T > &output_data, bool verbose=false) const
 Resample a volume into an existing output volume. More...
 
void resample (const carto::Volume< T > &input_data, const aims::AffineTransformation3d &transform, const T &background, const Point3df &output_location, T &output_value, int timestep) const
 Resample a volume at a single location. More...
 
void resample_inv (const carto::Volume< T > &input_data, const soma::Transformation3d &inverse_transform_to_mm, const T &background, const Point3df &output_location, T &output_value, int timestep) const
 Resample a volume at a single location. More...
 
void resample_inv (const carto::Volume< T > &input_data, const soma::Transformation3d &inverse_transform_to_mm, const T &background, carto::Volume< T > &output_data, bool verbose=false) const
 Resample a volume into an existing output volume. More...
 
void resample_inv_to_vox (const carto::Volume< T > &input_data, const soma::Transformation3d &inverse_transform_to_vox, const T &background, const Point3df &output_location, T &output_value, int timestep) const
 Resample a volume at a single location. More...
 
void setRef (const carto::rc_ptr< carto::Volume< T > > &ref)
 Set the input data to be resampled by the doit() methods. More...
 
const carto::rc_ptr< carto::Volume< T > > & ref () const
 Input data to be resampled by the doit() methods. More...
 
void setDefaultValue (T val)
 Set the background value to be used by the doit() methods. More...
 
const T & defaultValue () const
 Background value used by the doit() methods. More...
 
- Public Member Functions inherited from carto::RCObject
 RCObject ()
 
 RCObject (const RCObject &)
 
RCObjectoperator= (const RCObject &)
 
virtual ~RCObject ()
 

Protected Member Functions

void doResample (const carto::Volume< T > &inVolume, const soma::Transformation3d &transform3d, const T &outBackground, const Point3df &outLocation, T &outValue, int t) const CARTO_OVERRIDE
 Resample a volume at a single location. More...
 
virtual void doResampleChannel (const carto::Volume< ChannelType > &inVolume, const soma::Transformation3d &transform3d, const ChannelType &outBackground, const Point3df &outLocation, ChannelType &outValue, int t) const
 
void updateParameters (const carto::Volume< T > &inVolume, int t, bool verbose) const CARTO_OVERRIDE
 Update the cache of spline coefficients if needed. More...
 
virtual void updateParametersChannel (const carto::Volume< ChannelType > &inVolume, int t, bool verbose) const
 
void iirConvolveMirror (std::vector< double > &data) const
 
int getFold (int i, int size) const
 This method returns a mirror index when needed. More...
 
virtual double getBSplineWeight (int i, double x) const =0
 Returns $ B^n( x - i )$. More...
 

Protected Attributes

std::vector< double > _poles
 
double _gain
 
carto::VolumeRef< double > _splineCoefficients
 
const carto::Volume< ChannelType > * _lastvolume
 
int _lasttime
 
- Protected Attributes inherited from aims::Resampler< T >
carto::rc_ptr< carto::Volume< T > > _ref
 
_defval
 

Friends

template<bool , typename >
struct MultiChannelResamplerSwitch
 

Detailed Description

template<class T>
class aims::SplineResampler< T >

B-Spline-based resampling.

This is the base class for all resamplers based on B-spline interpolation as described by Unser, Aldroubi & Eden in IEEE PAMI (1991). A class computing the actual spline coefficient is derived for each spline order (1 to 7).

The resampling API is described in the base class, Resampler.

See also
LinearResampler CubicResampler QuadraticResampler QuinticResampler SixthOrderResampler SeventhOrderResampler

Coefficients are computed through a recursive scheme from an input reference volume. This recursive algorithm is fast for whole volume resampling, but slower for single points. A cache mechanism has thus been implemented so that coefficients are not recomputed if the input volume did not change.

For an input volume inVolume of dimensions [dimX, dimY, dimZ], coefficients are double precision values stored in a volume of same dimensions. Any value $ \mathrm{inVolume}(x, y, z)$ can then be interpolated through the formula:

\[ \mathrm{inVolume}(x, y, z) = \sum_i ( c_i * B^n(x - x_i) * B^n(y - y_i) * B^n(z - z_i) ) \]

with $ c_i$ the coefficients, $ B^n$ a $ n$-th order B-spline and $ (x_i, y_i, z_i)$ the discrete input image points. Given that the support of $ B^n$ is $ ]-(n+1)/2, (n+1)/2[$, only the $ (n+1)^d$ (with $ d$ the volume dimension in $[1,2,3]$) closest coefficients / discrete points are used to interpolate a given point.

Template Parameters
TVoxel data type

Definition at line 88 of file splineresampler.h.

Member Typedef Documentation

◆ ChannelType

Definition at line 139 of file splineresampler.h.

Constructor & Destructor Documentation

◆ SplineResampler()

template<class T >
aims::SplineResampler< T >::SplineResampler ( )

Definition at line 56 of file splineresampler_d.h.

◆ ~SplineResampler()

Member Function Documentation

◆ doResample()

template<class T >
void aims::SplineResampler< T >::doResample ( const carto::Volume< T > &  input_data,
const soma::Transformation3d inverse_transform,
const T &  background,
const Point3df output_location,
T &  output_value,
int  timestep 
) const
protectedvirtual

Resample a volume at a single location.

Parameters
[in]input_datadata to be resampled (its voxel size is not* taken into account)
[in]inverse_transformtransformation from output coordinates to coordinates of the input volume (unit: voxels)
[in]backgroundvalue set if the transformed point is outside of the input volume
[in]output_locationcoordinates in output space (source space of transform)
[out]output_valuevariable to be filled with resampled data
[in]timestepfor 4D volume, time step to be used

Implements aims::Resampler< T >.

Definition at line 768 of file splineresampler_d.h.

Referenced by aims::SplineResampler< T >::resample_inv_to_vox().

◆ doResampleChannel()

◆ getBSplineWeight()

template<class T>
virtual double aims::SplineResampler< T >::getBSplineWeight ( int  i,
double  x 
) const
protectedpure virtual

◆ getFold()

template<class T >
int aims::SplineResampler< T >::getFold ( int  i,
int  size 
) const
protected

This method returns a mirror index when needed.

Only inVolume's size spline coefficient are computed since "outside" coefficients are equal to their mirror inside the image domain. This method computes this mirror correspondance.

  • If i is in [0, size-1]: returns i
  • If i < 0: returns -1
  • If i >= size: returns size - (i - size) - 2

Definition at line 523 of file splineresampler_d.h.

References aims::mirrorCoeff(), and aims::SplineResampler< T >::resample_channel_inv_to_vox().

Referenced by aims::LinearResampler< T >::doResampleChannel(), aims::MajorityLabelResampler< T >::doResampleChannel(), aims::MedianResampler< T >::doResampleChannel(), and aims::SplineResampler< T >::doResampleChannel().

◆ getOrder()

◆ getSplineCoef()

template<class T >
carto::VolumeRef< double > aims::SplineResampler< T >::getSplineCoef ( const carto::Volume< T > &  inVolume,
int  t = 0,
bool  verbose = false 
)

Computes spline coefficients corresponding to an input volume.

Spline coefficients are recomputed only if one of these conditions is satisfied:

  • inVolume is different from the last volume used for coefficients computation (in the sense that they share the same adress)
  • t is different from the last t used for coefficients computation
  • A call to reset() has previously been made
Parameters
inVolumeinput image
tvolume to use in the T dimension in the case where inVolume is a time series.
verboseprint progress on stdout
Returns
Volume of double containing the coefficients in the image domain. Border coefficient need to be retrieved by mirror.

This method actually calls updateParameters() and returns the coeff container

Definition at line 72 of file splineresampler_d.h.

References aims::SplineResampler< T >::_splineCoefficients, and aims::SplineResampler< T >::updateParameters().

Referenced by aims::SplineResampler< T >::~SplineResampler().

◆ iirConvolveMirror()

template<class T >
void aims::SplineResampler< T >::iirConvolveMirror ( std::vector< double > &  data) const
protected

◆ resample_channel_inv_to_vox()

◆ resample_inv_to_vox()

template<class T >
void aims::SplineResampler< T >::resample_inv_to_vox ( const carto::Volume< T > &  input_data,
const soma::Transformation3d inverse_transform_to_vox,
const T &  background,
carto::Volume< T > &  output_data,
bool  verbose = false 
) const
virtual

Resample a volume into an existing output volume.

Parameters
[in]input_datadata to be resampled (its voxel size is not** taken into account)
[in]inverse_transform_to_voxtransformation from coordinates of the output volume (unit: mm), to coordinates of the input volume (unit: voxel)**
[in]backgroundvalue set in output regions that are outside of the transformed input volume
[in,out]output_dataexisting volume to be filled with resampled data (its pre-existing dimensions and voxel size are used)
[in]verboseprint progress to stdout

The transformations, referentials, and referential header attributes of output_data are not touched; it is up to the calling code to update them accordingly.

This method does not use the instance state set by setRef() or setDefaultValue().

Derived classes can override this method to optimize interpolation of a full volume. The base class method simply calls doResample for each point.

Reimplemented from aims::Resampler< T >.

Definition at line 751 of file splineresampler_d.h.

References aims::SplineResampler< T >::doResample().

Referenced by aims::SplineResampler< T >::updateParameters(), and aims::SplineResampler< T >::~SplineResampler().

◆ reset()

template<class T >
void aims::SplineResampler< T >::reset ( )

Clear the cache.

After a call to this method, a call to getSplineCoef() recomputes the spline coefficients in any case.

Definition at line 63 of file splineresampler_d.h.

References aims::SplineResampler< T >::_lasttime, and aims::SplineResampler< T >::_lastvolume.

Referenced by aims::SplineResampler< T >::~SplineResampler().

◆ updateParameters()

template<class T >
void aims::SplineResampler< T >::updateParameters ( const carto::Volume< T > &  inVolume,
int  t,
bool  verbose 
) const
protectedvirtual

Update the cache of spline coefficients if needed.

This method is called by all the resampling methods (resample() and doit()) before they call doResample for a given time point t.

Reimplemented from aims::Resampler< T >.

Definition at line 739 of file splineresampler_d.h.

References aims::SplineResampler< T >::resample_inv_to_vox().

Referenced by aims::SplineResampler< T >::getSplineCoef().

◆ updateParametersChannel()

Friends And Related Function Documentation

◆ MultiChannelResamplerSwitch

template<class T>
template<bool , typename >
friend struct MultiChannelResamplerSwitch
friend

Definition at line 137 of file splineresampler.h.

Member Data Documentation

◆ _gain

◆ _lasttime

◆ _lastvolume

template<class T>
const carto::Volume<ChannelType>* aims::SplineResampler< T >::_lastvolume
mutableprotected

◆ _poles

◆ _splineCoefficients


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