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

Resampling of data from a volume, applying a transformation. More...

#include <aims/resampling/resampler.h>

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

Public Member Functions

 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...
 
virtual 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
 Resample a volume into an existing output volume. 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

virtual void 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 =0
 Resample a volume at a single location. More...
 
virtual void updateParameters (const carto::Volume< T > &, int, bool) const
 Method called before doResample() for each time step. More...
 

Protected Attributes

carto::rc_ptr< carto::Volume< T > > _ref
 
_defval
 

Additional Inherited Members

- Public Types inherited from carto::RCObject
typedef int RefCounterType
 

Detailed Description

template<class T>
class aims::Resampler< T >

Resampling of data from a volume, applying a transformation.

The doit() and resample() methods can be used to apply an affine transformation (aims::AffineTransformation3d). They take a direct transformation, i.e. the transformation goes from the space of the input image (unit: mm) to the space of the output image (unit: mm). The transformation is inverted and normalized internally as needed, because the resamplers "pull" data by transforming output coordinates into input coordinates.

The doit() methods work on input data passed to the setRef() method. setDefaultValue() can also be called to set the background value.

The resample() methods provide stateless alternatives.

You can also use arbitrary non-affine transformations (inheriting soma::Transformation3d) by using the resample_inv() family of methods. In this case, you must pass the backward transformation (from output space to input space), because of the "pulling" mechanism described above.

Beware that contrary to the other methods, the resample_inv_to_vox() overloads take a transformation that maps to voxel coordinates of the input image. These methods can be slightly faster than resample_inv() because they map directly to the API that the actual resamplers are implementing (doResample()). The performance gain is especially noticeable for repeated calls to the overloads that perform resampling for a single point only.

Definition at line 75 of file resampler.h.

Constructor & Destructor Documentation

◆ Resampler()

template<typename T >
aims::Resampler< T >::Resampler ( )

Definition at line 61 of file resampler_d.h.

◆ ~Resampler()

template<class T>
virtual aims::Resampler< T >::~Resampler ( )
inlinevirtual

Definition at line 80 of file resampler.h.

Member Function Documentation

◆ defaultValue()

template<class T>
const T& aims::Resampler< T >::defaultValue ( ) const
inline

Background value used by the doit() methods.

Definition at line 318 of file resampler.h.

◆ doit() [1/2]

template<typename T>
void aims::Resampler< T >::doit ( const aims::AffineTransformation3d transform,
carto::Volume< T > &  output_data 
) const

Resample the input volume set with setRef() into an existing volume.

Parameters
[in]transformtransformation from coordinates of the input volume (unit: mm), to coordinates of the output* volume (unit: mm) (its inverse is used for resampling)
[in,out]output_dataexisting volume to be filled with resampled data (its pre-existing dimensions and voxel size are used)
Exceptions
std::runtime_errorif no input volume has been set with setRef.

The background value (to be used for regions that are outside of the input volume) can be set with setDefaultValue().

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

The level of verbosity is taken from carto::verbose (i.e. the --verbose command-line argument is honoured).

Definition at line 222 of file resampler_d.h.

Referenced by BlockMatching< T >::doit(), and aims::Resampler< float >::~Resampler().

◆ doit() [2/2]

template<typename T>
carto::VolumeRef< T > aims::Resampler< 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.

Parameters
[in]transformtransformation from coordinates of the input* volume (unit: mm), to coordinates of the output volume (unit: mm) (its inverse is used for resampling)
[in]dimX,dimY,dimZdimensions of the newly allocated volume
[in]voxel_sizevoxel size of the newly allocated volume (unit: mm)
Returns
a newly allocated volume containing the resampled data (its size along the t axis is the same as the input volume).
Exceptions
std::runtime_errorif no input volume has been set with setRef.

The background value (to be used for regions that are outside of the input volume) can be set with setDefaultValue().

The level of verbosity is taken from carto::verbose (i.e. the --verbose command-line argument is honoured).

The transformations, referentials, and referential header attributes of the new volume are reset if transform.isIdentity() is false:

  • the referential attribute is removed
  • each transformation in transformations is composed with transform so that the output volume still points to the original space. If that is not possible (e.g. the transformations attribute is missing or invalid), then a new transformation is added that points to the input volume.

Definition at line 233 of file resampler_d.h.

◆ doResample()

template<class T>
virtual void aims::Resampler< 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
protectedpure virtual

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

Implemented in aims::SplineResampler< T >, and aims::NearestNeighborResampler< T >.

Referenced by aims::Resampler< float >::defaultValue(), and aims::Resampler< float >::resample_inv_to_vox().

◆ ref()

template<class T>
const carto::rc_ptr<carto::Volume<T> >& aims::Resampler< T >::ref ( ) const
inline

Input data to be resampled by the doit() methods.

Definition at line 312 of file resampler.h.

Referenced by aims::Resampler< float >::resample_inv_to_vox().

◆ resample() [1/2]

template<typename T>
void aims::Resampler< T >::resample ( const carto::Volume< T > &  input_data,
const aims::AffineTransformation3d transform,
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 taken into account)
[in]transformtransformation from coordinates of the input volume (unit: mm), to coordinates of the output* volume (unit: mm) (its inverse is used for resampling)
[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.

Definition at line 175 of file resampler_d.h.

Referenced by aims::Resampler< float >::~Resampler().

◆ resample() [2/2]

template<typename T>
void aims::Resampler< T >::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.

Parameters
[in]input_datadata to be resampled (its voxel size is taken into account)
[in]transformtransformation from coordinates of the input* volume (unit: mm), to output coordinates (its inverse is used for resampling)
[in]backgroundvalue set in output regions that are outside of the transformed input volume
[in]output_locationcoordinates in output space (destination space of transform)
[out]output_valuevariable to be filled with resampled data
[in]timestepfor 4D volume, time step to be used

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

Definition at line 198 of file resampler_d.h.

◆ resample_inv() [1/2]

template<typename T>
void aims::Resampler< T >::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.

Parameters
[in]input_datadata to be resampled (its voxel size is taken into account)
[in]inverse_transformtransformation from output coordinates to coordinates of the input volume (unit: mm)
[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

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

Definition at line 151 of file resampler_d.h.

Referenced by aims::Resampler< float >::~Resampler().

◆ resample_inv() [2/2]

template<typename T>
void aims::Resampler< T >::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.

Parameters
[in]input_datadata to be resampled (its voxel size is taken into account)
[in]transformtransformation from coordinates of the output* volume (unit: mm), to coordinates of the input volume *(unit: mm)*
[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().

Definition at line 129 of file resampler_d.h.

◆ resample_inv_to_vox() [1/2]

template<class T>
void aims::Resampler< T >::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
inline

Resample a volume at a single location.

Parameters
[in]input_datadata to be resampled (its voxel size is not** taken into account)
[in]inverse_transform_to_voxtransformation from output coordinates to coordinates of the input volume (unit: voxel)
[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

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

Definition at line 266 of file resampler.h.

Referenced by aims::Resampler< float >::resample_inv_to_vox().

◆ resample_inv_to_vox() [2/2]

template<typename T>
void aims::Resampler< 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 in aims::SplineResampler< T >.

Definition at line 69 of file resampler_d.h.

◆ setDefaultValue()

template<class T>
void aims::Resampler< T >::setDefaultValue ( val)
inline

Set the background value to be used by the doit() methods.

Definition at line 315 of file resampler.h.

◆ setRef()

template<typename T>
void aims::Resampler< T >::setRef ( const carto::rc_ptr< carto::Volume< T > > &  ref)

Set the input data to be resampled by the doit() methods.

Definition at line 307 of file resampler_d.h.

Referenced by BlockMatching< T >::doit(), and aims::Resampler< float >::resample_inv_to_vox().

◆ updateParameters()

template<class T>
virtual void aims::Resampler< T >::updateParameters ( const carto::Volume< T > &  ,
int  ,
bool   
) const
inlineprotectedvirtual

Method called before doResample() for each time step.

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

The base class version of this method does nothing.

Reimplemented in aims::SplineResampler< T >.

Definition at line 349 of file resampler.h.

Referenced by aims::Resampler< float >::resample_inv_to_vox().

Member Data Documentation

◆ _defval

template<class T>
T aims::Resampler< T >::_defval
protected

◆ _ref

template<class T>
carto::rc_ptr<carto::Volume<T> > aims::Resampler< T >::_ref
protected

Definition at line 351 of file resampler.h.

Referenced by aims::Resampler< float >::ref().


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