NPL
Neurological Programs and Libraries
MRImage Functions

Functions

ptr< MRImage > npl::createMRImage (size_t ndim, const size_t *size, PixelT type)
 Creates a new MRImage with dimensions set by ndim, and size set by size. Output pixel type is decided by type variable. More...
 
ptr< MRImage > npl::createMRImage (const std::vector< size_t > &dim, PixelT type)
 Creates a new MRImage with dimensions set by ndim, and size set by size. Output pixel type is decided by type variable. More...
 
ptr< MRImage > npl::createMRImage (size_t ndim, const size_t *size, PixelT type, void *ptr, std::function< void(void *)> deleter)
 Creates a new MRImage with dimensions set by ndim, and size set by size. Output pixel type is decided by type variable. More...
 
ptr< MRImage > npl::createMRImage (const std::vector< size_t > &dim, PixelT type, void *ptr, std::function< void(void *)> deleter)
 Creates a new MRImage with dimensions set by ndim, and size set by size. Output pixel type is decided by type variable. More...
 
std::ostream & npl::operator<< (std::ostream &out, const MRImage &img)
 Writes out information about an MRImage. More...
 
void npl::gaussianSmooth1D (ptr< MRImage > inout, size_t dim, double stddev)
 Gaussian smooths an image in 1 direction. Questionable whether it works. Seems to shift image. More...
 
ptr< MRImage > npl::smoothDownsample (ptr< const MRImage > in, double sigma, double spacing=-1)
 Performs smoothing in each dimension, then downsamples so that pixel spacing is roughly equal to FWHM. More...
 
ptr< MRImage > npl::resample (ptr< const MRImage > in, double *spacing, double(*window)(double, double)=hannWindow)
 Performs fourier resampling using fourier transform and the provided window function. More...
 
ptr< MRImage > npl::shiftImage (ptr< MRImage > in, size_t len, double *vect)
 Uses fourier shift theorem to shift an image. More...
 
void npl::writeComplex (std::string basename, ptr< const MRImage > in, bool absPhase=false)
 Writes a pair of images, one real, one imaginary or if absPhase is set to true then an absolute image and a phase image. More...
 
ptr< MRImage > npl::fft_forward (ptr< const MRImage > in, const std::vector< size_t > &in_osize)
 Performs forward FFT transform in N dimensions. More...
 
ptr< MRImage > npl::fft_backward (ptr< const MRImage > in, const std::vector< size_t > &in_osize)
 Performs inverse FFT transform in N dimensions. More...
 
int npl::rotateImageShearKern (ptr< MRImage > inout, double rx, double ry, double rz, double(*kern)(double, double)=npl::lanczosKern)
 Rotates an image around the center using shear decomposition followed by kernel-based shearing. Rotation order is Rz, Ry, Rx, and about the center of the image. This means that 1D interpolation will be used. More...
 
int npl::rotateImageShearFFT (ptr< MRImage > inout, double rx, double ry, double rz, double(*window)(double, double)=npl::sincWindow)
 Rotates an image around the center using shear decomposition followed by FFT-based shearing. Rotation order is Rz, Ry, Rx, and about the center of the image. More...
 
ptr< MRImage > npl::rigidTransform (ptr< MRImage > in, double rx, double ry, double rz, double sx, double sy, double sz)
 Rigid Transforms an image. More...
 
ptr< MRImage > npl::diffOfGauss (ptr< const MRImage > in, double sd1, double sd2)
 Computes difference of gaussians. More...
 
double npl::overlapRatio (ptr< const MRImage > a, ptr< const MRImage > b)
 Computes the overlap of the two images' in 3-space. More...
 
ptr< MRImage > npl::resampleNN (ptr< const MRImage > in, ptr< const MRImage > atlas, PixelT type=UNKNOWN_TYPE)
 Performs nearest neighbor resasmpling of input to atlas. More...
 
ptr< MRImage > npl::resampleNN (ptr< const MRImage > in, double *newspace, PixelT type=UNKNOWN_TYPE)
 Performs nearest neighbor resasmpling of input to atlas. More...
 
ptr< MRImage > npl::randImage (PixelT type, double mean, double sd, size_t x, size_t y, size_t z, size_t t)
 Create random image, with gaussian distribution. More...
 

Detailed Description

Function Documentation

ptr<MRImage> npl::createMRImage ( size_t  ndim,
const size_t *  size,
PixelT  type 
)

Creates a new MRImage with dimensions set by ndim, and size set by size. Output pixel type is decided by type variable.

Parameters
ndimnumber of image dimensions
sizesize of image, in each dimension
typePixel type npl::PixelT
Returns
New image, default orientation
ptr<MRImage> npl::createMRImage ( const std::vector< size_t > &  dim,
PixelT  type 
)

Creates a new MRImage with dimensions set by ndim, and size set by size. Output pixel type is decided by type variable.

Parameters
dimsize of image, in each dimension, number of dimensions decied by length of size vector
typePixel type npl::PixelT
Returns
New image, default orientation
ptr<MRImage> npl::createMRImage ( size_t  ndim,
const size_t *  size,
PixelT  type,
void *  ptr,
std::function< void(void *)>  deleter 
)

Creates a new MRImage with dimensions set by ndim, and size set by size. Output pixel type is decided by type variable.

Parameters
ndimnumber of image dimensions
sizesize of image, in each dimension
typePixel type npl::PixelT
ptrPointer to data block
deleterfunction to delete data block
Returns
New image, default orientation
ptr<MRImage> npl::createMRImage ( const std::vector< size_t > &  dim,
PixelT  type,
void *  ptr,
std::function< void(void *)>  deleter 
)

Creates a new MRImage with dimensions set by ndim, and size set by size. Output pixel type is decided by type variable.

Parameters
dimsize of image, in each dimension, number of dimensions decied by length of size vector
typePixel type npl::PixelT
ptrPointer to data block
deleterfunction to delete data block
Returns
New image, default orientation
ptr<MRImage> npl::diffOfGauss ( ptr< const MRImage in,
double  sd1,
double  sd2 
)

Computes difference of gaussians.

Parameters
inGaussian smooths an image twice and subtracts
sd1Standard deviation of first image
sd2Standard deviation of second image
Returns
Difference of two different gaussians
ptr<MRImage> npl::fft_backward ( ptr< const MRImage in,
const std::vector< size_t > &  in_osize 
)

Performs inverse FFT transform in N dimensions.

Parameters
inInput image
in_osizeSize of output image. If this is smaller than the input then the frequency domain will be trunkated, if it is larger then the fourier domain will be padded ( output upsampled )
Returns
Frequency domain of input. Note the output will be COMPLEX128/CDOUBLE type
ptr<MRImage> npl::fft_forward ( ptr< const MRImage in,
const std::vector< size_t > &  in_osize 
)

Performs forward FFT transform in N dimensions.

Parameters
inInput image
in_osizeSize of output image (will be padded up to this prior to FFT)
Returns
Frequency domain of input. Note the output will be COMPLEX128/CDOUBLE type
void npl::gaussianSmooth1D ( ptr< MRImage inout,
size_t  dim,
double  stddev 
)

Gaussian smooths an image in 1 direction. Questionable whether it works. Seems to shift image.

Parameters
inoutInput/Output image
dimDirection to smooth in
stddevin real space, for example millimeters.
std::ostream & npl::operator<< ( std::ostream &  out,
const MRImage img 
)

Writes out information about an MRImage.

Parameters
outOutput ostream
imgImage to write information about
Returns
More ostream
double npl::overlapRatio ( ptr< const MRImage a,
ptr< const MRImage b 
)

Computes the overlap of the two images' in 3-space.

Parameters
aImage
bImage
Returns
Ratio of b that overlaps with a's grid
ptr<MRImage> npl::randImage ( PixelT  type,
double  mean,
double  sd,
size_t  x,
size_t  y,
size_t  z,
size_t  t 
)

Create random image, with gaussian distribution.

Parameters
typeType of pixels
meanmean of Gaussian
sdStandard deviation of distribution
xX size
yY size
zZ size
tTime size, default 0 time
Returns
Output MRImage
ptr<MRImage> npl::resample ( ptr< const MRImage in,
double *  spacing,
double(*)(double, double)  window = hannWindow 
)

Performs fourier resampling using fourier transform and the provided window function.

Parameters
inInput image
spacingDesired output spacing
windowWindow function to reduce ringing
Returns
Smoothed and downsampled image
ptr<MRImage> npl::resampleNN ( ptr< const MRImage in,
ptr< const MRImage atlas,
PixelT  type = UNKNOWN_TYPE 
)

Performs nearest neighbor resasmpling of input to atlas.

Parameters
inInput image
atlas
typeof output pixel (defaults to input type)
Returns
Input image resampled into atlas space
ptr<MRImage> npl::resampleNN ( ptr< const MRImage in,
double *  newspace,
PixelT  type = UNKNOWN_TYPE 
)

Performs nearest neighbor resasmpling of input to atlas.

Parameters
inInput image
newspacenew spacing
typepixel type of output (defaults to input)
Returns
Input image resampled into atlas space
ptr<MRImage> npl::rigidTransform ( ptr< MRImage in,
double  rx,
double  ry,
double  rz,
double  sx,
double  sy,
double  sz 
)

Rigid Transforms an image.

Parameters
inInput/output image
rxRotation about x axis
ryRotation about y axis
rzRotation about z axis
sxshift about x axis
syshift about y axis
szshift about z axis
int npl::rotateImageShearFFT ( ptr< MRImage inout,
double  rx,
double  ry,
double  rz,
double(*)(double, double)  window = npl::sincWindow 
)

Rotates an image around the center using shear decomposition followed by FFT-based shearing. Rotation order is Rz, Ry, Rx, and about the center of the image.

Parameters
inoutInput/output image
rxRotation about x axis
ryRotation about y axis
rzRotation about z axis
rzRotation about z axis
windowWindow function to apply in fourier domain
int npl::rotateImageShearKern ( ptr< MRImage inout,
double  rx,
double  ry,
double  rz,
double(*)(double, double)  kern = npl::lanczosKern 
)

Rotates an image around the center using shear decomposition followed by kernel-based shearing. Rotation order is Rz, Ry, Rx, and about the center of the image. This means that 1D interpolation will be used.

Parameters
inoutInput/output image
rxRotation about x axis
ryRotation about y axis
rzRotation about z axis
kernKernel to use for interplation
ptr<MRImage> npl::shiftImage ( ptr< MRImage in,
size_t  len,
double *  vect 
)

Uses fourier shift theorem to shift an image.

Parameters
inInput image to shift
lenlength of dx array
vectmovement in physical coordinates, will be rotated using image orientation prior to shifting
Returns
shifted image
ptr<MRImage> npl::smoothDownsample ( ptr< const MRImage in,
double  sigma,
double  spacing = -1 
)

Performs smoothing in each dimension, then downsamples so that pixel spacing is roughly equal to FWHM.

Parameters
inInput image
sigmaStandard deviation for smoothing
spacingOuptut image spacing (isotropic). If this is <= 0, then sigma will be used, which is a very conservative downsampling.
Returns
Smoothed and downsampled image
void npl::writeComplex ( std::string  basename,
ptr< const MRImage in,
bool  absPhase = false 
)

Writes a pair of images, one real, one imaginary or if absPhase is set to true then an absolute image and a phase image.

Parameters
basenameBase filename _abs.nii.gz and _phase.nii.gz or _re.nii.gz and _im.nii.gz will be appended, depending on absPhase
inInput image
absPhaseWhether the break up into absolute and phase rather than re/imaginary