Methods list

Methods usage

SPECTrecon.PlanPSFType
PlanPSF{T,Tf,Ti}( ; nx::Int, nz::Int, px::Int, pz::Int, T::Type)

Struct for storing work arrays and factors for 2D convolution for one thread. Each PSF is px × pz

  • T datatype of work arrays (subtype of AbstractFloat)
  • nx::Int = 128 (ny implicitly the same)
  • nz::Int = nx image size is [nx,nx,nz]
  • px::Int = 1
  • pz::Int = px (PSF size)
  • padsize::Tuple : (padup, paddown, padleft, padright)
  • workmat [nx+padup+paddown, nz+padleft+padright] 2D padded image for FFT convolution
  • workvecx [nx+padup+paddown,]: 1D work vector
  • workvecz [nz+padleft+padright,]: 1D work vector
  • img_compl [nx+padup+paddown, nz+padleft+padright]: 2D [complex] padded image for FFT
  • ker_compl [nx+padup+paddown, nz+padleft+padright]: 2D [complex] padded image for FFT
  • fft_plan::Tf plan for doing FFT; see plan_fft!
  • ifft_plan::Ti plan for doing IFFT; see plan_ifft!
source
SPECTrecon.PlanRotateType
PlanRotate{T, R}

Struct for storing work arrays and factors for 2D square image rotation for one thread

  • T datatype of work arrays (default Float32)
  • R::RotateMode
  • nx::Int image size
  • padsize::Int : pad each side of image by this much
  • workmat1 [nx + 2 * padsize, nx + 2 * padsize] padded work matrix
  • workmat2 [nx + 2 * padsize, nx + 2 * padsize] padded work matrix
  • interp::SparseInterpolator{T, 2, 1}
source
SPECTrecon.RotateModeType
RotateMode

Type to control image rotation method.

  • Rotate1D : 3-pass 1D interpolation
  • Rotate2D : 1-pass 2D bilinear interpolation
source
SPECTrecon.SPECTplanType
SPECTplan

Struct for storing key factors for a SPECT system model

  • T datatype of work arrays
  • imgsize size of image
  • px,pz psf dimension
  • imgr [nx, ny, nz] 3D rotated version of image
  • add_img [nx, ny, nz] 3D image for adding views and backprojection
  • mumap [nx,ny,nz] attenuation map, must be 3D, possibly zeros()
  • mumapr [nx, ny, nz] 3D rotated mumap
  • exp_mumapr [nx, nz] 2D exponential rotated mumap
  • psfs [px,pz,ny,nview] point spread function, must be 4D, with px andpz` odd, and symmetric for each slice
  • nview number of views, must be integer
  • viewangle set of view angles, must be from 0 to 2π
  • interpmeth interpolation method: :one means 1d; :two means 2d
  • mode pre-allocation method: :fast means faster; :mem means use less memory
  • dy voxel size in y direction (dx is the same value)
  • nthread number of CPU threads used to process data; must be integer
  • planrot Vector of struct PlanRotate
  • planpsf Vector of struct PlanPSF

Currently code assumes the following:

  • each of the nview projection views is [nx,nz]
  • nx = ny
  • uniform angular sampling
  • psf is symmetric
  • multiprocessing using # of threads specified by Threads.nthreads()
source
Base.showMethod
show(io::IO, mime::MIME"text/plain", vp::Vector{<:PlanPSF})
source
Base.showMethod
show(io::IO, mime::MIME"text/plain", vp::Vector{<:PlanRotate})
source
Base.showMethod
show(io::IO, ::MIME"text/plain", plan::PlanRotate)
source
Base.showMethod
show(io::IO, ::MIME"text/plain", plan::SPECTplan)
source
Base.sizeofMethod
sizeof(::SPECTplan)

Show size in bytes of SPECTplan object.

source
Base.sizeofMethod
sizeof(::PlanRotate)

Show size in bytes of PlanRotate object.

source
SPECTrecon.AblockMethod
Ablock(plan, nblocks)

Generate a vector of linear maps for OSEM. -plan: SPECTrecon plan -nblocks: Number of blocks in OSEM

source
SPECTrecon.backproject!Method
backproject!(image, views, plan ; index)

Backproject multiple views into image. Array image is not initialized to zero; caller must do that.

source
SPECTrecon.backprojectMethod
image = backproject(views, mumap, psfs, dy; interpmeth, kwargs...)

SPECT backproject views using attenuation map mumap and PSF array psfs for pixel size dy. This method initializes the plan as a convenience. Most users should use backproject! instead after initializing those, for better efficiency.

source
SPECTrecon.backprojectMethod
image = backproject(views, plan ; kwargs...)

SPECT backproject views; this allocates the returned 3D array.

source
SPECTrecon.fft_conv!Method
fft_conv!(output, image3, ker3, plans)

Mutating version of convolving a 3D image3 with a set of 2D symmetric kernels stored in 3D array ker3 using foreach.

source
SPECTrecon.fft_conv!Method
fft_conv!(output, img, ker, plan)

Convolve 2D image img with 2D (symmetric!) kernel ker using FFT, storing the result in output.

source
SPECTrecon.fft_conv_adj2!Method
fft_conv_adj2!(output, image2, ker3, plans)

Mutating version of adjoint of convolving a 2D image2 with each 2D kernel in the 3D array ker3.

source
SPECTrecon.imfilterz!Method
imfilterz!(plan)

FFT-based convolution of plan.img_compl and kernel plan.ker_compl (not centered), storing result in plan.workmat.

source
SPECTrecon.imrotate!Function
imrotate!(output, image3, θ, plans, ntasks=nthreads)

In-place version of rotating a 3D image3 by θ in counter-clockwise direction (opposite to imrotate in Julia) using foreach with ntasks.

source
SPECTrecon.imrotate!Method
imrotate!(output, image3, θ, plans, :thread)

Mutating version of rotating a 3D image3 by θ in counter-clockwise direction (opposite of imrotate in Julia) writing result into output, using Threads.@threads.

source
SPECTrecon.imrotate!Method
imrotate!(output, img, θ, plan)

Rotate a 2D image img by angle θ ∈ [0,2π] in counter-clockwise direction (opposite to imrotate in Julia) using 3-pass 1d linear interpolations.

source
SPECTrecon.imrotate!Method
imrotate!(output, img, θ, plan)

Rotate a 2D image img by angle θ ∈ [0,2π] in counter-clockwise direction (opposite to imrotate in Julia) using 2d bilinear interpolation.

source
SPECTrecon.imrotateMethod
imrotate(img, θ; method::Symbol=:two)

Rotate a 2D image img by angle θ ∈ [0,2π] in counter-clockwise direction (opposite to imrotate in Julia) using either 2d linear interpolation (for :two) or 3-pass 1D interpolation (for :one)

source
SPECTrecon.imrotate_adj!Method
imrotate_adj!(output, img, θ, plan)

Adjoint of imrotate! for a 2D image using 3-pass 1d linear interpolations.

source
SPECTrecon.imrotate_adjMethod
imrotate_adj(img, θ; method::Symbol=:two)

Adjoint of rotating a 2D image img by angle θ ∈ [0,2π] in counter-clockwise direction (opposite to imrotate in Julia) using either 2d linear interpolations or 3-pass 1D interpolation.

source
SPECTrecon.linearinterp!Method
linearinterp!(A, s, e, len)

Assign key values in SparseInterpolator (linear) A that are calculated from s, e and len. s means start (x[1]) e means end (x[end]) len means length (length(x))

source
SPECTrecon.mlem!Method
mlem!(out, x0, ynoisy, background, A; niter = 20)

Inplace version of ML-EM algorithm for emission tomography image reconstruction.

  • out: Output
  • x0: Initial guess
  • ynoisy: (Noisy) measurements
  • background: Background effects, e.g., scatters
  • A: System matrix
  • niter: Number of iterations
source
SPECTrecon.mlemMethod
mlem(x0, ynoisy, background, A; niter = 20)

ML-EM algorithm for emission tomography image reconstruction.

  • x0: Initial guess
  • ynoisy: (Noisy) measurements
  • background: Background effects, e.g., scatters
  • A: System matrix
  • niter: Number of iterations
source
SPECTrecon.osem!Method
osem!(out, x0, ynoisy, background, Ab; niter = 20)

OS-EM algorithm for SPECT reconstruction.

  • out: Output
  • x0: Initial estimate
  • ynoisy: (Noisy) measurements
  • background: Background effects, e.g., scatters
  • Ab: Vector of system matrix
  • niter: Number of iterations
source
SPECTrecon.osemMethod
osem(x0, ynoisy, background, Ab; niter = 20)

OS-EM algorithm for SPECT reconstruction.

  • x0: Initial guess
  • ynoisy: (Noisy) measurements
  • background: Background effects, e.g., scatters
  • Ab: Vector of system matrix
  • niter: Number of iterations
source
SPECTrecon.pad2sizezero!Method
pad2sizezero!(output, img, padsize)

Non-allocating version of padding: `output[paddims[1]+1 : paddims[1]+dims[1], paddims[2]+1 : paddims[2]+dims[2]]) .= img

source
SPECTrecon.padzero!Method
padzero!(output, img, pad_x, pad_y)

Mutating version of padding a 2D image by filling zeros. Output has size (size(img, 1) + padsize[1] + padsize[2], size(img, 2) + padsize[3] + padsize[4]).

source
SPECTrecon.plan_psfMethod
plan_psf( ; nx::Int, nz::Int, px::Int, pz::Int, nthread::Int, T::Type)

Make Vector of structs for storing work arrays and factors for 2D convolution with SPECT depth-dependent PSF model, threaded across planes parallel to detector. Option

  • nx::Int = 128
  • nz::Int = nx
  • px::Int = 1
  • pz::Int = px PSF size is px × pz
  • T : datatype of work arrays, defaults to Float32
  • nthread::Int # of threads, defaults to Threads.nthreads()
source
SPECTrecon.plan_rotateMethod
plan_rotate(nx::Int; T::Type{<:AbstractFloat}, method::Symbol)

Make Vector of PlanRotate structs for storing work arrays and factors for threaded rotation of a stack of 2D square images.

Input

  • nx::Int must equal to ny (square images only)

Option

  • T : datatype of work arrays, defaults to Float32
  • method::Symbol : default is :two for 2D interpolation; use :one for 3-pass rotation with 1D interpolation
  • nthread::Int : default to Threads.nthreads() The other useful option is 1 when rotating just one image.
source
SPECTrecon.project!Method
project!(views, image, plan; index)

Project image into multiple views with indexes index (defaults to 1:nview). The 3D views array must be pre-allocated, but need not be initialized.

source
SPECTrecon.project!Method
project!(view, plan, image, thid, viewidx)

SPECT projection of image into a single view with index viewidx. The view must be pre-allocated but need not be initialized to zero.

source
SPECTrecon.project!Method
project!(view, plan, image, viewidx)

SPECT projection of image into a single view with index viewidx. The view must be pre-allocated but need not be initialized to zero.

source
SPECTrecon.projectMethod
views = project(image, mumap, psfs, dy; interpmeth, kwargs...)

Convenience method for SPECT forward projector that does all allocation including initializing plan.

In

  • image : 3D array (nx,ny,nz)
  • mumap : (nx,ny,nz) 3D attenuation map, possibly zeros()
  • psfs : 4D PSF array
  • dy::RealU : pixel size

Option

  • interpmeth : :one or :two
source
SPECTrecon.projectMethod
views = project(image, plan ; kwargs...)

Convenience method for SPECT forward projector that allocates and returns views.

source
SPECTrecon.psf_gaussMethod
psf_gauss( ; ny, px, pz, fwhm_start, fwhm_end, fwhm, fwhm_x, fwhm_z, T)

Create depth-dependent Gaussian PSFs having specified full-width half-maximum (FHWM) values.

Options

  • 'ny::Int = 128'
  • 'px::Int = 11' (should be odd)
  • 'pz::Int = px' (should be odd)
  • 'fwhm_start::Real = 1'
  • 'fwhm_end::Real = 4'
  • 'fwhm::AbstractVector{<:Real} = range(fwhmstart, fwhmend, ny)'
  • 'fwhm_x::AbstractVector{<:Real} = fwhm,
  • 'fwhmz::AbstractVector{<:Real} = fwhmx'
  • 'T::Type == Float32'

Returned psf is [px, pz, ny] where each PSF sums to 1.

source
SPECTrecon.rot180!Method
rot180!(B::AbstractMatrix, A::AbstractMatrix)

In place version of rot180, returning rotation of A in B.

source
SPECTrecon.rotl90!Method
rotl90!(B::AbstractMatrix, A::AbstractMatrix)

In place version of rotl90, returning rotation of A in B.

source
SPECTrecon.rotr90!Method
rotr90!(B::AbstractMatrix, A::AbstractMatrix)

In place version of rotr90, returning rotation of A in B.

source