ODIN
Classes | Enumerations | Functions
Classes of the ODIN data processing framework (odindata library)

Classes

class  ComplexData< N_rank >
 
class  Converter
 
struct  correlationResult
 
struct  fmriResult
 
class  Data< T, N_rank >
 
class  FileIO
 
class  FileFormat
 
struct  FileReadOpts
 
struct  FileWriteOpts
 
class  FilterChain
 
struct  fitpar
 
class  ModelFunction
 
class  FunctionFitInterface
 
class  FunctionFitDerivative
 
struct  ExponentialFunction
 
struct  ExponentialFunctionWithOffset
 
struct  GaussianFunction
 
struct  SinusFunction
 
struct  GammaVariateFunction
 
struct  PolynomialFunction< N_rank >
 
struct  LinearFunction
 
class  DownhillSimplex
 
class  FunctionFitDownhillSimplex
 
struct  GriddingPoint< N_rank >
 
class  Gridding< T, N_rank >
 
class  CoordTransformation< T, N_rank, OnPixelRot >
 
class  Image
 
class  ImageSet
 
class  Integrand
 
class  FunctionIntegral
 
struct  statisticResult
 
class  Step< T >
 
class  StepFactory< T >
 

Enumerations

enum  dataDim
 

Functions

template<typename T , int N_rank>
correlationResult correlation (const Array< T, N_rank > &x, const Array< T, N_rank > &y)
 
template<typename T , int N_rank>
correlationResult kendall (const Array< T, N_rank > &x, const Array< T, N_rank > &y)
 
fmriResult fmri_eval (const Data< float, 1 > &timecourse, const Data< float, 1 > &designvec)
 
template<int N_rank>
Array< float, N_rank > polyniomial_fit (const Array< float, N_rank > &value_map, const Array< float, N_rank > &reliability_map, unsigned int polynom_order, float kernel_size, bool only_zero_reliability=false)
 
Data< float, 1 > solve_linear (const Data< float, 2 > &A, const Data< float, 1 > &b, float sv_truncation=0.0)
 
ComplexData< 1 > solve_linear (const ComplexData< 2 > &A, const ComplexData< 1 > &b, float sv_truncation=0.0)
 
Data< float, 1 > eigenvalues (const Data< float, 2 > &A)
 
template<typename T , int N_rank>
statisticResult statistics (const Array< T, N_rank > &ensemble, const Array< T, N_rank > *mask=0)
 
template<typename T , int N_rank>
median (const Array< T, N_rank > &ensemble, const Array< T, N_rank > *mask=0)
 
template<typename T , int N_rank>
weightmean (const Array< T, N_rank > &ensemble, const Array< float, N_rank > &weight)
 
Data< float, 1 > unwrap_phase (const Data< float, 1 > &phase, int startindex=0)
 
template<typename T >
Array< T, 1 > matrix_product (const Array< T, 2 > &matrix, const Array< T, 1 > &vector)
 
template<typename T >
Array< T, 2 > matrix_product (const Array< T, 2 > &matrix1, const Array< T, 2 > &matrix2)
 
template<typename T >
Array< T, 1 > vector_product (const Array< T, 1 > &u, const Array< T, 1 > &v)
 
template<typename T , int N_rank>
bool operator== (const TinyVector< T, N_rank > &t1, const TinyVector< T, N_rank > &t2)
 
template<typename T , int N_rank>
bool operator!= (const TinyVector< T, N_rank > &t1, const TinyVector< T, N_rank > &t2)
 
template<typename T , int N_rows, int N_columns>
TinyVector< T, N_rows > operator* (const TinyMatrix< T, N_rows, N_columns > &matrix, const TinyVector< T, N_columns > &vector)
 
template<typename T , int N_rank>
bool same_shape (const Array< T, N_rank > &a1, const Array< T, N_rank > &a2, const TinyVector< int, N_rank > &dimmask=1)
 
template<int N_rank>
bool on_grid (const TinyVector< int, 3 > &shape, const TinyVector< int, 3 > &index)
 
template<typename T >
bool check_range (T &val, T min, T max)
 
template<typename T , int N_rank>
void clip_max (Data< T, N_rank > &data, T val)
 
template<typename T , int N_rank>
void clip_min (Data< T, N_rank > &data, T val)
 
bool PolynomialFunction< N_rank >::fit (const Array< float, 1 > &yvals, const Array< float, 1 > &ysigma, const Array< float, 1 > &xvals)
 
Array< float, 1 > PolynomialFunction< N_rank >::get_function (const Array< float, 1 > &xvals) const
 
Array< float, N_rank > Gridding< T, N_rank >::init (const TinyVector< int, N_rank > &dst_shape, const TinyVector< float, N_rank > &dst_extent, const STD_vector< GriddingPoint< N_rank > > &src_coords, const LDRfilter &kernel, float kernel_diameter)
 
template<int N_rank_in>
Array< T, N_rank > Gridding< T, N_rank >::operator() (const Array< T, N_rank_in > &src, unsigned int offset=0) const
 
 CoordTransformation< T, N_rank, OnPixelRot >::CoordTransformation (const TinyVector< int, N_rank > &shape, const TinyMatrix< float, N_rank, N_rank > &rotation, const TinyVector< float, N_rank > &offset, float kernel_size=2.5)
 
Array< T, N_rank > CoordTransformation< T, N_rank, OnPixelRot >::operator() (const Array< T, N_rank > &A) const
 

Detailed Description

ODIN data processing framework (odindata library) Input/Output of medical image data (odindata library)

Enumeration Type Documentation

◆ dataDim

enum dataDim

The 4 dimensions of data arrays used in FileIO and Filter components, these dimensions are

  • timeDim: Time direction
  • sliceDim: Slice direction
  • phaseDim: Phase direction
  • readDim: Read direction

Definition at line 45 of file fileio.h.

Function Documentation

◆ check_range()

template<typename T >
bool check_range ( T &  val,
min,
max 
)

Limits value 'val' to range (min,max)

Definition at line 197 of file utils.h.

◆ clip_max()

template<typename T , int N_rank>
void clip_max ( Data< T, N_rank > &  data,
val 
)

Clip all values above maximum value

Definition at line 211 of file utils.h.

◆ clip_min()

template<typename T , int N_rank>
void clip_min ( Data< T, N_rank > &  data,
val 
)

Clip all values below minimum value

Definition at line 223 of file utils.h.

◆ CoordTransformation()

template<typename T , int N_rank, bool OnPixelRot>
CoordTransformation< T, N_rank, OnPixelRot >::CoordTransformation ( const TinyVector< int, N_rank > &  shape,
const TinyMatrix< float, N_rank, N_rank > &  rotation,
const TinyVector< float, N_rank > &  offset,
float  kernel_size = 2.5 
)

Initializes a transformation of an array with shape 'shape' to new coordinates using a gridding algorithm. New coordinates (...z',y',x') are obtained by multiplying (...z,y,x) with 'rotation' and adding 'offset'. Origin is the center of the data set.

Definition at line 251 of file gridding.h.

◆ correlation()

template<typename T , int N_rank>
correlationResult correlation ( const Array< T, N_rank > &  x,
const Array< T, N_rank > &  y 
)

Linear correlation between vectors 'x' and 'y'

Definition at line 61 of file correlation.h.

◆ eigenvalues()

Data<float,1> eigenvalues ( const Data< float, 2 > &  A)

Computes the eigenvalues of symmetric matrix A using LAPACK (or with GSL if LAPACK is not available).

◆ fit()

template<int N_rank>
bool PolynomialFunction< N_rank >::fit ( const Array< float, 1 > &  yvals,
const Array< float, 1 > &  ysigma,
const Array< float, 1 > &  xvals 
)

polynomial fitting routine. Fits the function to the y-values 'yvals', and optionally the corresponding error bars 'ysigma' and x-values 'xvals'. If no error-bars are given they are all set to 1.0 and if no x-vals are given equidistant points with an increment of one are chosen, i.e. xvals(i)=i; Returns true on success.

Definition at line 362 of file fitting.h.

◆ fmri_eval()

fmriResult fmri_eval ( const Data< float, 1 > &  timecourse,
const Data< float, 1 > &  designvec 
)

Returns an fMRI analysis of 'timecourse' using 'designvec'

◆ get_function()

template<int N_rank>
Array< float, 1 > PolynomialFunction< N_rank >::get_function ( const Array< float, 1 > &  xvals) const

Returns the polynomial function values for x-values 'xvals' using the current polynomial coefficients.

Definition at line 401 of file fitting.h.

◆ init()

template<typename T , int N_rank>
Array< float, N_rank > Gridding< T, N_rank >::init ( const TinyVector< int, N_rank > &  dst_shape,
const TinyVector< float, N_rank > &  dst_extent,
const STD_vector< GriddingPoint< N_rank > > &  src_coords,
const LDRfilter kernel,
float  kernel_diameter 
)

Initializes a gridding object to perform a gridding operation with the following parameters:

  • dst_shape: The dimensions of the gridded array
  • dst_extent: The total extent of the gridded array, the gridded array is created symmetrically about the origin
  • src_coords: The coordinates of the input array: First is tthe coordinate, second is an extra weight for the coordinate
  • kernel: The gridding kernel used
  • kernel_diameter: The maximum diameter of the gridding kernel in units of the source coordinates Returns the density of source points on the gridded array

Definition at line 96 of file gridding.h.

◆ kendall()

template<typename T , int N_rank>
correlationResult kendall ( const Array< T, N_rank > &  x,
const Array< T, N_rank > &  y 
)

Kendall's rank-correlation between vectors 'x' and 'y' re-implemented from NRC, section 14.6

Definition at line 106 of file correlation.h.

◆ matrix_product() [1/2]

template<typename T >
Array<T,1> matrix_product ( const Array< T, 2 > &  matrix,
const Array< T, 1 > &  vector 
)

matrix-vector product

Definition at line 42 of file utils.h.

◆ matrix_product() [2/2]

template<typename T >
Array<T,2> matrix_product ( const Array< T, 2 > &  matrix1,
const Array< T, 2 > &  matrix2 
)

matrix-matrix product

Definition at line 72 of file utils.h.

◆ median()

template<typename T , int N_rank>
T median ( const Array< T, N_rank > &  ensemble,
const Array< T, N_rank > *  mask = 0 
)

Returns the median of 'ensemble'.

Definition at line 146 of file statistics.h.

◆ on_grid()

template<int N_rank>
bool on_grid ( const TinyVector< int, 3 > &  shape,
const TinyVector< int, 3 > &  index 
)

Checks whether multidimensional index 'index' is within the array with shape 'shape'

Definition at line 181 of file utils.h.

◆ operator!=()

template<typename T , int N_rank>
bool operator!= ( const TinyVector< T, N_rank > &  t1,
const TinyVector< T, N_rank > &  t2 
)

Unequal-comparison operator for TinyVectors

Definition at line 139 of file utils.h.

◆ operator()() [1/2]

template<typename T , int N_rank, bool OnPixelRot>
Array< T, N_rank > CoordTransformation< T, N_rank, OnPixelRot >::operator() ( const Array< T, N_rank > &  A) const

Transforms 'A' to new coordinates and returns the result

Definition at line 283 of file gridding.h.

◆ operator()() [2/2]

template<typename T , int N_rank>
template<int N_rank_in>
Array< T, N_rank > Gridding< T, N_rank >::operator() ( const Array< T, N_rank_in > &  src,
unsigned int  offset = 0 
) const

Put input array 'src' on the grid by stepping linearly through the inidices. Gridding will start at the linear index 'offset' of the 'src_coords' specified in the init() function. Returns gridded array.

Definition at line 194 of file gridding.h.

◆ operator*()

template<typename T , int N_rows, int N_columns>
TinyVector<T,N_rows> operator* ( const TinyMatrix< T, N_rows, N_columns > &  matrix,
const TinyVector< T, N_columns > &  vector 
)

Product operator for TinyMatrix*TinyMatrix (defined here because the product(TinyMatrix,TinyMatrix) function is missing in blitz-0.10)

Definition at line 149 of file utils.h.

◆ operator==()

template<typename T , int N_rank>
bool operator== ( const TinyVector< T, N_rank > &  t1,
const TinyVector< T, N_rank > &  t2 
)

Equal-comparison operator for TinyVectors

Definition at line 129 of file utils.h.

◆ polyniomial_fit()

template<int N_rank>
Array<float,N_rank> polyniomial_fit ( const Array< float, N_rank > &  value_map,
const Array< float, N_rank > &  reliability_map,
unsigned int  polynom_order,
float  kernel_size,
bool  only_zero_reliability = false 
)

Fits an N_rank-dimensional polynomial of order 'polynom_order' to each point of the array using the values of its neighbours regarding their reliability (i.e. their relative weight for the fit). Parameters are:

  • value_map: The array to be fitted
  • reliability_map: The reliability of each point
  • polynom_order: Order of the polynom
  • kernel_size: Size of the neighbourhood of the pixel which is considered for the fit (using a Gaussian kernel with this FWHM)
  • only_zero_reliability: Fit only pixel with zero reliabiliy

This function returns the fitted array

Definition at line 554 of file fitting.h.

◆ same_shape()

template<typename T , int N_rank>
bool same_shape ( const Array< T, N_rank > &  a1,
const Array< T, N_rank > &  a2,
const TinyVector< int, N_rank > &  dimmask = 1 
)

Compares array shapes while discarding dimensions with zero value in dimmask. Returns true if equal.

Definition at line 168 of file utils.h.

◆ solve_linear() [1/2]

ComplexData<1> solve_linear ( const ComplexData< 2 > &  A,
const ComplexData< 1 > &  b,
float  sv_truncation = 0.0 
)

Solves the complex linear system A x = b and returns x. All singular values less than sv_truncation * max_singular_value will be set to zero. The algorithm uses the singular value decomposition from LAPACK (or from GSL if LAPACK is not available).

◆ solve_linear() [2/2]

Data<float,1> solve_linear ( const Data< float, 2 > &  A,
const Data< float, 1 > &  b,
float  sv_truncation = 0.0 
)

Solves the linear system A x = b and returns x. All singular values less than sv_truncation * max_singular_value will be set to zero. The algorithm uses the singular value decomposition from LAPACK (or from GSL if LAPACK is not available).

◆ statistics()

template<typename T , int N_rank>
statisticResult statistics ( const Array< T, N_rank > &  ensemble,
const Array< T, N_rank > *  mask = 0 
)

Returns a statistical description of 'ensemble'. If mask is non-zero, only values with mask!=0 are considered.

Definition at line 79 of file statistics.h.

◆ unwrap_phase()

Data<float,1> unwrap_phase ( const Data< float, 1 > &  phase,
int  startindex = 0 
)

unwrap phase in one dimension, starting at position 'startindex'

◆ vector_product()

template<typename T >
Array<T,1> vector_product ( const Array< T, 1 > &  u,
const Array< T, 1 > &  v 
)

cross (vector) product (always in 3 dimensions)

Definition at line 109 of file utils.h.

◆ weightmean()

template<typename T , int N_rank>
T weightmean ( const Array< T, N_rank > &  ensemble,
const Array< float, N_rank > &  weight 
)

Returns the weighted mean of 'ensemble' using 'weight' for the weighting.

Definition at line 180 of file statistics.h.