![]() |
aimsalgo
5.1.2
Neuroimaging image processing
|
B-Spline-based resampling. More...
#include <aims/resampling/splineresampler.h>


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 | setVerboseStream (std::ostream &stream) |
| write verbose progress information to this stream. More... | |
| 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 &) | |
| RCObject & | operator= (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 | |
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 |
| T | _defval |
| std::ostream * | _verbose_stream |
Friends | |
| template<bool , typename > | |
| struct | MultiChannelResamplerSwitch |
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.
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
can then be interpolated through the formula:
with
the coefficients,
a
-th order B-spline and
the discrete input image points. Given that the support of
is
, only the
(with
the volume dimension in
) closest coefficients / discrete points are used to interpolate a given point.
| T | Voxel data type |
Definition at line 88 of file splineresampler.h.
| typedef carto::DataTypeTraits<T>::ChannelType aims::SplineResampler< T >::ChannelType |
Definition at line 139 of file splineresampler.h.
| aims::SplineResampler< T >::SplineResampler |
Definition at line 56 of file splineresampler_d.h.
|
inline |
Definition at line 93 of file splineresampler.h.
|
protectedvirtual |
Resample a volume at a single location.
| [in] | input_data | data to be resampled (its voxel size is not* taken into account) |
| [in] | inverse_transform | transformation from output coordinates to coordinates of the input volume (unit: voxels) |
| [in] | background | value set if the transformed point is outside of the input volume |
| [in] | output_location | coordinates in output space (source space of transform) |
| [out] | output_value | variable to be filled with resampled data |
| [in] | timestep | for 4D volume, time step to be used |
Implements aims::Resampler< T >.
Definition at line 781 of file splineresampler_d.h.
|
protectedvirtual |
Reimplemented in aims::MedianResampler< T >, aims::MajorityLabelResampler< T >, and aims::LinearResampler< T >.
Definition at line 84 of file splineresampler_d.h.
References carto::RawConverter< class, class >::convert(), carto::VolumeProxy< class >::getSizeX(), carto::VolumeProxy< class >::getSizeY(), carto::VolumeProxy< class >::getSizeZ(), std::numeric_limits< class >::max(), and soma::Transformation3d::transform().
Referenced by aims::MultiChannelResamplerSwitch< false, T >::doResample(), and aims::MultiChannelResamplerSwitch< true, T >::doResample().
|
protectedpure virtual |
Returns
.
This method is defined by each nth order derived class
Implemented in aims::SixthOrderResampler< T >, aims::SeventhOrderResampler< T >, aims::QuinticResampler< T >, aims::QuarticResampler< T >, aims::QuadraticResampler< T >, aims::MedianResampler< T >, aims::MajorityLabelResampler< T >, aims::LinearResampler< T >, and aims::CubicResampler< T >.
|
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.
i is in [0, size-1]: returns ii < 0: returns -1i >= size: returns size - (i - size) - 2 Definition at line 537 of file splineresampler_d.h.
References aims::mirrorCoeff().
|
pure virtual |
Spline order (1 to 7)
Implemented in aims::SixthOrderResampler< T >, aims::SeventhOrderResampler< T >, aims::QuinticResampler< T >, aims::QuarticResampler< T >, aims::QuadraticResampler< T >, aims::MedianResampler< T >, aims::MajorityLabelResampler< T >, aims::LinearResampler< T >, and aims::CubicResampler< 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| inVolume | input image |
| t | volume to use in the T dimension in the case where inVolume is a time series. |
| verbose | print progress on stdout |
This method actually calls updateParameters() and returns the coeff container
Definition at line 72 of file splineresampler_d.h.
References verbose.
|
protected |
Definition at line 421 of file splineresampler_d.h.
References EPSILON.
| void aims::SplineResampler< T >::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 |
Definition at line 547 of file splineresampler_d.h.
References carto::Volume< class >::at(), carto::VolumeProxy< class >::getSizeT(), carto::VolumeProxy< class >::getSizeX(), carto::VolumeProxy< class >::getSizeY(), carto::VolumeProxy< class >::getSizeZ(), carto::VolumeProxy< class >::getVoxelSize(), aims::ProgressInfo< class, class >::print(), and verbose.
Referenced by aims::MultiChannelResamplerSwitch< true, T >::resample_inv_to_vox().
|
virtual |
Resample a volume into an existing output volume.
| [in] | input_data | data to be resampled (its voxel size is not** taken into account) |
| [in] | inverse_transform_to_vox | transformation from coordinates of the output volume (unit: mm), to coordinates of the input volume (unit: voxel)** |
| [in] | background | value set in output regions that are outside of the transformed input volume |
| [in,out] | output_data | existing volume to be filled with resampled data (its pre-existing dimensions and voxel size are used) |
| [in] | verbose | print 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 764 of file splineresampler_d.h.
References verbose.
| 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.
|
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 753 of file splineresampler_d.h.
References verbose.
|
protectedvirtual |
Definition at line 260 of file splineresampler_d.h.
References carto::Volume< class >::at(), carto::Converter< class, class >::convert(), carto::VolumeProxy< class >::getSize(), carto::VolumeProxy< class >::getVoxelSize(), carto::ref< class >::release(), carto::rc_ptr< class >::reset(), and verbose.
Referenced by aims::MultiChannelResamplerSwitch< true, T >::doResample(), and aims::MultiChannelResamplerSwitch< false, T >::updateParameters().
Definition at line 137 of file splineresampler.h.
|
protected |
Definition at line 190 of file splineresampler.h.
Referenced by aims::CubicResampler< T >::CubicResampler(), aims::QuadraticResampler< T >::QuadraticResampler(), aims::QuarticResampler< T >::QuarticResampler(), aims::QuinticResampler< T >::QuinticResampler(), aims::SeventhOrderResampler< T >::SeventhOrderResampler(), and aims::SixthOrderResampler< T >::SixthOrderResampler().
|
mutableprotected |
Definition at line 195 of file splineresampler.h.
|
mutableprotected |
Definition at line 194 of file splineresampler.h.
|
protected |
Definition at line 189 of file splineresampler.h.
Referenced by aims::CubicResampler< T >::CubicResampler(), aims::QuadraticResampler< T >::QuadraticResampler(), aims::QuarticResampler< T >::QuarticResampler(), aims::QuinticResampler< T >::QuinticResampler(), aims::SeventhOrderResampler< T >::SeventhOrderResampler(), and aims::SixthOrderResampler< T >::SixthOrderResampler().
|
mutableprotected |
Definition at line 193 of file splineresampler.h.