Methods list

Methods usage

Sinograms.CtFanArcMethod
CtFanArc(::Val{:ge1} ; kwargs...)

GE Lightspeed system CT geometry.

option

  • unit::RealU = 1 or use 1mm
  • (see CtFanArc)

out

  • CtFanArc

These numbers are published in IEEE T-MI Oct. 2006, p.1272-1283 wang:06:pwl.

julia> CtFanArc(Val(:ge1))
CtFanArc{Float64, Float32, CtSourceCircle} :
 ns::Int64 888
 nt::Int64 64
 ds::Float64 1.0239
 dt::Float64 1.0964
 offset_s::Float32 1.25
 offset_t::Float32 0.0
 na::Int64 984
 orbit::Float32 360.0
 orbit_start::Float32 0.0
 source_offset::Float64 0.0
 dsd::Float64 949.075
 dod::Float64 408.075
 src::CtSourceCircle CtSourceCircle()
source
Sinograms.CtFanArcMethod
CtFanArc( ; ns nt ds dt offset_s offset_t
    na orbit orbit_start
    dsd = 4ns * ds, dod = ns * ds)
CtFanArc(:short ; ...)

Constructor with named keywords. See ?CtGeom for documentation.

  • Use :short argument to specify a short scan, in which case na will be scaled down proportionally as well.
julia> CtFanArc()
CtFanArc{Float32, Float32, CtSourceCircle} :
 ns::Int64 128
 nt::Int64 64
 ds::Float32 1.0
 dt::Float32 1.0
 offset_s::Float32 0.0
 offset_t::Float32 0.0
 na::Int64 64
 orbit::Float32 360.0
 orbit_start::Float32 0.0
 source_offset::Float32 0.0
 dsd::Float32 512.0
 dod::Float32 128.0
 src::CtSourceCircle CtSourceCircle()
source
Sinograms.CtFanFlatMethod
CtFanFlat( ; ns nt ds dt offset_s offset_t
    na orbit orbit_start
    dsd = 4ns * ds, dod = ns * ds)
CtFanFlat(:short ; ...)

Constructor with named keywords. See ?CtGeom for documentation.

  • Use :short argument to specify a short scan, in which case na will be scaled down proportionally as well.
julia> CtFanFlat()
CtFanFlat{Float32, Float32, CtSourceCircle} :
 ns::Int64 128
 nt::Int64 64
 ds::Float32 1.0
 dt::Float32 1.0
 offset_s::Float32 0.0
 offset_t::Float32 0.0
 na::Int64 64
 orbit::Float32 360.0
 orbit_start::Float32 0.0
 source_offset::Float32 0.0
 dsd::Float32 512.0
 dod::Float32 128.0
 src::CtSourceCircle CtSourceCircle()
source
Sinograms.CtGeomType
CtGeom{Td,To,Ts}

Abstract type for representing ray geometries for 3D CT imaging.

The projection view coordinates are (s,t) where s denotes the transaxial sampling and t denotes the axial direction (along z).

Common fields

  • ns size of each projection view
  • nt
  • ds detector pixel spacing
  • dt
  • offset_s unitless detector center offset (usually 0 or 0.25)
  • offset_t unitless, usually 0
  • na # of projection views, aka nϕ or nβ
  • orbit source orbit in degrees (or Unitful)
  • orbit_start starting angle in degrees (or Unitful) orbit and orbit_start must both be unitless (degrees) or have same units.
  • src::CtSource describes the X-ray CT source trajectory. Primary support for CtSourceCircle().

Additional fields for CtFan types:

  • dsd distance from source to detector
  • dod distance from origin to detector
  • source_offset usually 0

Units:

  • ds, dt, source_offset, dsd, dod must all be unitless or have the same units.

Basic methods

  • angles (na) in degrees
  • dims (ns, nt, na)
  • ones = ones(Float32, ns,nt,na)
  • zeros = zeros(Float32, ns,nt,na)
  • rays iterator of (u,v,ϕ,θ) samples
  • downsample(st, down) reduce sampling by integer factor
  • oversample(st, over)
  • ct_geom_plot3 plot system geometry

Non-exported helper functions for developers:

  • _s (ns) s sample locations
  • _t (nt) t sample locations
  • _ws = (ns-1)/2 + offset_s "middle" sample position
  • _wt = (nt-1)/2 + offset_t
  • _ar (na) source angles [radians]
  • _rfov max radius within FOV
  • _zfov axial FOV
  • _xds (nb) center of detector elements (beta=0)
  • _yds (nb) ""
  • _tau(rg, x, y) projected s/ds for each (x,y) pair (length(x), na)
  • _shape(rg, proj [,:]) reshape proj into array (ns,nt,na[,:])
  • _unitv(rg [, (is,it,ia)]) unit 'vector' with single nonzero element

For fan beam:

  • _dso = dsd - dod distance from source to origin (Inf for parallel beam)
  • _dfs distance from source to detector focal spot (0 for 3rd gen CT, Inf for flat detectors)
  • _gamma(rg [,s]) (ns) gamma sample values radians, optionally given s values
  • _gamma_max = max(|γ|) half of fan angle radians, if offset_s == 0
  • _cone_angle (half) cone angle on axis (s=0): +/- angle

Notes

Use sino_geom() instead for 2D geometries.

source
Sinograms.CtParMethod
CtPar( ; ns nt ds dt offset_s offset_t na orbit orbit_start)

Constructor with named keywords. See ?CtGeom for documentation.

julia> CtPar()
CtPar{Float32, Float32, CtSourceCircle} :
 ns::Int64 128
 nt::Int64 64
 ds::Float32 1.0
 dt::Float32 1.0
 offset_s::Float32 0.0
 offset_t::Float32 0.0
 na::Int64 60
 orbit::Float32 180.0
 orbit_start::Float32 0.0
 src::CtSourceCircle CtSourceCircle()
source
Sinograms.CtSourceHelixType
CtSourceHelix

Type for representing helical X-ray CT source trajectories having constant pitch.

Fields

  • pitch helix pitch (unitless fraction of zFOV, default 0)
  • source_z0 initial z position of x-ray source (default 0) It should have same units as detector pixels etc.
source
Sinograms.FBPNormalPlanType
FBPNormalPlan{S,I,H,P}

Struct type for storing "normal" FBP plan.

The view_weight can include (products of)

  • Parker weighting for short scans
  • view-wise weighting from fbp_sino_weight todo
  • weighting for possibly nonuniform angles from _angle_weights
source
Sinograms.FDKplanType
FDKplan{C,I,H,V}

Struct type for storing FDK plan.

The view_weight can include (products of)

  • Parker weighting for short scans
  • view-wise CBCT weighting from fdk_weight_cyl
  • weighting for possibly nonuniform angles from _angle_weights
source
Sinograms.SinoFanArcMethod
SinoFanArc(Val(:ge1) ; kwargs...)

GE Lightspeed system CT geometry.

option

  • unit::RealU = 1 or use 1mm
  • nb::Int = 888 # of detector channels
  • d::RealU = 1.0239 channel spacing
  • offset::Real = 1.25 for "quarter-detector" offset
  • na::Int = 984 # of angles
  • orbit::Union{Symbol,Real} = 360 use :short for short scan
  • dsd::RealU = 949.075
  • dod::RealU = 408.075

out

  • SinoFanArc

These numbers are published in IEEE T-MI Oct. 2006, p.1272-1283 wang:06:pwl.

julia> SinoFanArc(Val(:ge1))
SinoFanArc{Float32, Float32} :
 nb::Int64 888
 d::Float32 1.0239
 offset::Float32 1.25
 na::Int64 984
 orbit::Float32 360.0
 orbit_start::Float32 0.0
 source_offset::Float32 0.0
 dsd::Float32 949.075
 dod::Float32 408.075
source
Sinograms.SinoFanArcMethod
SinoFanArc( ; nb d offset na orbit orbit_start
    source_offset, dsd = 4 * nb * d, dod = nb * d)
SinoFanArc(:short ; ...)

Constructor with named keywords. See ?SinoGeom for documentation.

  • d, source_offset, dsd, dod must all have the same units.

  • Use :short argument to specify a short scan, in which case na will be scaled down proportionally as well.

julia> SinoFanArc()
SinoFanArc{Float32, Float32} :
 nb::Int64 128
 d::Float32 1.0
 offset::Float32 0.0
 na::Int64 200
 orbit::Float32 360.0
 orbit_start::Float32 0.0
 source_offset::Float32 0.0
 dsd::Float32 512.0
 dod::Float32 128.0
source
Sinograms.SinoFanFlatMethod
SinoFanFlat( ; nb d offset na orbit orbit_start
    source_offset, dsd= 4 * nb * d, dod = nb * d)
SinoFanFlat(:short ; ...)

Constructor with named keywords. See ?SinoGeom for documentation.

  • Use :short argument to specify a short scan, in which case na will be scaled down proportionally as well.
julia> SinoFanFlat()
SinoFanFlat{Float32, Float32} :
 nb::Int64 128
 d::Float32 1.0
 offset::Float32 0.0
 na::Int64 200
 orbit::Float32 360.0
 orbit_start::Float32 0.0
 source_offset::Float32 0.0
 dsd::Float32 512.0
 dod::Float32 128.0
source
Sinograms.SinoGeomType
SinoGeom{Td,To}

Abstract type for representing ray geometries for 2D sinograms. This describes the sampling characteristics of a given sinogram for a 2D parallel or fan-beam system.

Common fields

  • nb # of "radial" samples, aka nr or ns
  • d aka dr or ds, "radial" sample spacing
  • offset unitless sample offset (usually 0 or 0.25)
  • na # of angular samples, aka or default: 2 * floor(Int, nb * π/2 / 2)
  • orbit source orbit in degrees (or Unitful)
  • orbit_start starting angle in degrees (or Unitful) orbit and orbit_start must both be unitless (degrees) or have same units.

Additional fields for SinoFan types

  • dsd distance from source to detector
  • dod distance from origin to detector
  • source_offset usually 0

Units:

  • d, source_offset, dsd, dod must all have the same units.

Basic methods

  • angles (na) in degrees
  • dims (nb, na)
  • ones = ones(Float32, nb, na)
  • zeros = zeros(Float32, nb, na)
  • rays iterator of (r, ϕ) parallel-beam coordinate tuples of size (nb, na)
  • downsample(st, down) reduce sampling by integer factor
  • oversample(st, over)
  • sino_geom_plot! plot system geometry

Non-exported helper functions for developers:

  • _ds|dr radial sample spacing (NaN for :moj)
  • _s (nb) s sample locations
  • _w = (nb-1)/2 + offset "middle" sample position
  • _ar (na) source angles [radians]
  • _rfov radial FOV
  • _xds (nb) center of detector elements (beta=0)
  • _yds (nb) ""
  • _tau(rg, x, y) projected s/ds for each (x,y) pair (length(x), na)
  • _shape(rg, sino [,:]) reshape sino into array (nb,na[,:])
  • _unitv(rg [, (ib,ia)]) unit 'vector' with single nonzero element

For mojette:

  • _d_ang (na) angle-dependent radial spacing

For fan beam:

  • _dso = dsd - dod distance from source to origin (Inf for parallel beam)
  • _dfs distance from source to detector focal spot (0 for 3rd gen CT, Inf for flat detectors)
  • _gamma(rg [,s]) (nb) gamma sample values radians, optionally given s values
  • _gamma_max = max(|γ|) half of fan angle radians, if offset_s == 0

Notes

  • Use ct_geom() instead for 3D axial or helical cone-beam CT.
source
Sinograms.SinoMojMethod
SinoMoj( ; nb d offset na orbit orbit_start)

Constructor with named keywords. See ?SinoGeom for documentation.

julia> SinoMoj()
SinoMoj{Float32, Float32} :
 nb::Int64 128
 d::Float32 1.0
 offset::Float32 0.0
 na::Int64 200
 orbit::Float32 180.0
 orbit_start::Float32 0.0
source
Sinograms.SinoParMethod
SinoPar( ; nb d offset na orbit orbit_start)

Constructor with named keywords. See ?SinoGeom for documentation.

julia> SinoPar()
SinoPar{Float32, Float32} :
 nb::Int64 128
 d::Float32 1.0
 offset::Float32 0.0
 na::Int64 200
 orbit::Float32 180.0
 orbit_start::Float32 0.0
source
Sinograms.WindowType
Window{S,T}

Data type for FBP apodizing windows, with given window shape::S and cutoff::T.

julia> Window(Hamming(), 0.8)
Window{Hamming, Float64}(Hamming(), 0.8)
source
Sinograms.WindowVectType
WindowVect{V} <: AbstractWindowShape

A user-specified window vector, constructed via WindowVect(v::V), where v is a AbstractVector. Caution: length(v) must be appropriate for the padded sinogram size.

source
Sinograms._angle_weightsMethod
_angle_weights(ar::AbstractVector{<:RealU})

Angular weighting to # scale projections by dβ (aka dϕ or da) for Riemann-like integration.

input

  • ar (na,) angles in radians.

output

For now, simply the scalar (ar[begin+1]-ar[begin]) / n180 where n180 is the number of full multiples of 180°. This simplifies to π/na for typical 180° and 360° scans. For a fan-beam short scan (180° + fan angle), it also simplifies to π/na, where the "excess" is handled by Parker weighting.

Examples

julia> rg = SinoPar(); Sinograms._angle_weights(Sinograms._ar(rg)), π/rg.na
(0.015707962f0, 0.015707963267948967)
julia> rg = SinoPar(;orbit=360); Sinograms._angle_weights(Sinograms._ar(rg)), π/rg.na
(0.015707962f0, 0.015707963267948967)
julia> rg = SinoFanArc(:short); Sinograms._angle_weights(Sinograms._ar(rg)), π/rg.na
(0.04842342f0, 0.04487989505128276)

For now, only equally spaced views are supported, but this is where unequal spacing would be handled, and the output would be of size (na,).

source
Sinograms._dfsMethod

_dfs(::Union{SinoFan,CtFan}) Distance from detector arc focal spot to source for fan-beam.

source
Sinograms._dsoMethod
_dso(rg::Union{SinoFan,CtFan}) = rg.dsd - rg.dod

Distance from source to origin for fan-beam.

source
Sinograms._gammaFunction
_gamma(rg::RayGeom [, s])

Return gamma (γ: fan angle) values for fan-beam geometry.

source
Sinograms._realeMethod
_reale(x ; rtol = …)

Return the real part of x, but warn if the imaginary part is too large. By default this uses a 100× larger rtol than isapprox due to FFT errors.

source
Sinograms._showMethod
_show(io::IO, ::MIME"text/plain", st::Any)

Informative way to show fields of a struct (composite type).

source
Sinograms._tauMethod
_tau(rg::RayGeom, x, y)

Projected s/ds, useful for footprint center and support. Returns Matrix of size length(x) × rg.na.

source
Sinograms._unitvFunction
_unitv([T=Float32], rg:RayGeom [, pos::Tuple])

Projection views with a single non-zero ray value at position pos (default: middle).

source
Sinograms._xdsMethod
_xds(rg::RayGeom)

Center x positions of detectors (for beta = 0), for central row of detector.

source
Sinograms._ydsMethod
_yds(rg::RayGeom)

Center y positions of detectors (for beta = 0), for central row of detector.

source
Sinograms.anglesMethod
angles(rg::RayGeom) =

Return vector of angles for this ray geometry, in whatever units the user used to specify orbit and orbit_start, typically degrees.

source
Sinograms.backproject_bddMethod
backproject_bdd(sinogram::AbstractMatrix{<:T}, geo)

Generates a reconstructed image using the back projection algorithm, for a sinogram and a tuple of geometry definitions

source
Sinograms.cb_arc_to_parMethod
u, v, azim, polar = cb_arc_to_par(s, t, β, dso, dod)

Convert from cone-beam arc (3rd gen CT) coordinates to parallel-beam coordinates.

in

  • s,t detector coordinates
  • β X-ray source angle, measured counter-clockwise from y axis
  • dso,dsd distances for the geometry

out

  • u,v transaxial and axial parallel-beam detector coordinates
  • azim ϕ transaxial or azimuthal angle (radians)
  • polar θ polar angle (radians)
source
Sinograms.cb_flat_to_parMethod
u, v, azim, polar = cb_flat_to_par(s, t, β, dso, dod)

Convert from cone-beam flat panel coordinages to parallel-beam coordinates.

in

  • s,t detector coordinates
  • β X-ray source angle, measured counter-clockwise from y axis
  • dso,dsd distances for the geometry

out

  • u,v transaxial and axial parallel-beam detector coordinates
  • azim ϕ transaxial or azimuthal angle (radians)
  • polar θ polar angle (radians)
source
Sinograms.cbct_backMethod
cbct_back(proj, rg, ig)

Cone-beam backprojector for feldkamp.jl

in

  • proj (ns,nt,na) cone-beam projection views
  • rg::CtGeom
  • ig::ImageGeom

out

  • img (nx,ny,nz) back-projection result
source
Sinograms.cbct_back_fan!Method
cbct_back_fan!(...)

This should work even for non-uniformly spaced source angles (for a circular source trajectory).

source
Sinograms.downsampleMethod
downsample(rg, down::Int)
downsample(rg, down::NTuple{3,Int})

Down-sample CT geometry (for testing with small problems).

source
Sinograms.fbpFunction
fbp(plan, sino)

Filtered back-projection (FBP) reconstruction.

in

  • plan::FBPplan
  • sino::AbstractArray{<:Number} (nb,na,...) sinogram(s)

out

  • image::Matrix{<:Number} (nx,ny,...) reconstructed image(s)
source
Sinograms.fbp!Method
fbp!(image, sino ; orbit::Real = 180, kwargs...)

FBP reconstruction from a parallel-beam sinogram of size [nr × nϕ]. Writes result into image matrix.

Input

  • sino::AbstractMatrix

Options for SinoPar constructor

  • dr : sinogram radial spacing; default 1
  • orbit : angular range in degrees; default 180
  • orbit_start : angular range in degrees; default 0

Options for ImageGeom

  • dx, dy, deltas, offset_x, offset_y, offsets
  • rmax maximum radius for mask

Options

  • kwargs : passed to plan_fbp

Output

  • image::AbstractMatrix is mutated
source
Sinograms.fbpMethod
fbp(sino::AbstractMatrix ; kwargs...)

FBP reconstruction from a parallel-beam sinogram of size [nr × nϕ]. Returns an image of size [nx × ny].

Options

  • nx : default nr
  • ny : default nx
  • kwargs : passed to fbp!
source
Sinograms.fbp_backFunction
img = fbp_back(rg, ig, sino ; ia_skip)

2D pixel-driven backprojection for FBP.

in

  • rg::SinoGeom
  • ig::ImageGeom
  • sino::AbstractArray{<:Number} sinogram(s) (line integrals), usually ramp filtered

option

  • ia_skip::Int downsample in angle to save time for quick tests (default: 1)

out

  • img::Array{<:Number} reconstructed image(s)
source
Sinograms.fbp_back_fan!Method
fbp_back_fan!(image, sino, sinβ, cosβ,
    dsd, dso, source_offset, is_arc,
    ds, offset, xc, yc, mask ; ia_skip)

Mutating version of pixel-driven back-projection for a grid of (xc,yc) pixel center locations for sinogram sino from a fan-beam geometry. It uses Threads. It assumes the angles are equally spaced over [0,π).

in

  • sino::Matrix{<:Number} (nb, na) usually ramp-filtered
  • sinβ::Vector{<:Real} (na)
  • cosβ::Vector{<:Real} (na)
  • dsd,dso,source_offset::RealU geometry
  • is_arc::Bool
  • ds::RealU ray spacing
  • offset::Real detector offset (usually 0)
  • xc::Vector{<:RealU} (nx) pixel centers
  • yc::Vector{<:RealU} (ny) pixel centers
  • mask::Matrix{Bool} (nx, ny) which pixels to reconstruct

option

  • ia_skip::Int default 1

out

  • image::Matrix (nx, ny) matrix to be mutated
source
Sinograms.fbp_back_fanMethod
fbp_back_fan(sino, betas,
    dsd, dso, dfs, source_offset, is_arc,
    ds, offset, xc, yc, mask ; ia_skip, T)

Pixel-driven back-projection for a grid of (xc,yc) pixel center locations for sinogram sino from a fan-beam geometry. It assumes the angles are equally spaced over [0,π).

in

  • sino::Matrix{<:Number} (nb, na) usually ramp-filtered
  • betas::Vector{<:Real} (na) in radians
  • dsd,dso,source_offset::RealU geometry
  • is_arc::Bool arc or flat?
  • ds::RealU ray spacing
  • offset::Real detector offset (usually 0)
  • xc::Vector{<:RealU} (nx) pixel centers
  • yc::Vector{<:RealU} (ny) pixel centers
  • mask::Matrix{Bool} (nx, ny) which pixels to reconstruct

option

  • ia_skip::Int default 1
  • T::Type{<:Number} usually same as eltype(sino)

out

  • image::Matrix (nx, ny)
source
Sinograms.fbp_back_fan_xyMethod
fbp_back_fan_xy(sino, sinβ, cosβ,
    dsd_ds, dso_ds, source_offset_ds, is_arc,
    wb, x, y ; T)

Pixel-driven back-projection for a single (x,y) location for sinogram sino from a fan-beam geometry. It assumes the angles are equally spaced over [0,π).

in

  • sino::Matrix{<:Number} (nb, na) usually ramp-filtered
  • sinβ::Vector{<:Real} (na)
  • cosβ::Vector{<:Real} (na)
  • dsd_ds,dso_ds,source_offset_ds::Real geometry, normalized
  • wb::Real = (nb+1)/2 + offset where usually offset=0
  • x,y::Real pixel center location, normalized by ray spacing

option

  • T::Type{<:Number} typically same as eltype(sino)

out

  • Returns a scalar of type T.
source
Sinograms.fbp_back_par!Method
fbp_back_par!(image, sino, sinϕ, cosϕ,
    ds, offset, xc, yc, mask ; ia_skip)

Mutating version of pixel-driven back-projection for a grid of (xc,yc) pixel center locations for sinogram sino from a parallel-beam geometry. It uses Threads. It assumes the angles are equally spaced over [0,π).

in

  • sino::Matrix{<:Number} (nb, na) usually ramp-filtered
  • sinϕ::Vector{<:Real} (na)
  • cosϕ::Vector{<:Real} (na)
  • ds::RealU ray spacing
  • offset::Real detector offset (usually 0)
  • xc::Vector{<:RealU} (nx) pixel centers
  • yc::Vector{<:RealU} (ny) pixel centers
  • mask::Matrix{Bool} (nx, ny) which pixels to reconstruct

option

  • ia_skip::Int default 1

out

  • image::Matrix (nx, ny) matrix to be mutated
source
Sinograms.fbp_back_parMethod
fbp_back_par(sino, angles,
    ds, offset, xc, yc, mask ; ia_skip, T)

Pixel-driven back-projection for a grid of (xc,yc) pixel center locations for sinogram sino from a parallel-beam geometry. It assumes the angles are equally spaced over [0,π).

in

  • sino::Matrix{<:Number} (nb, na) usually ramp-filtered
  • angles::Vector{<:Real} (na) in radians
  • ds::RealU ray spacing
  • offset::Real detector offset (usually 0)
  • xc::Vector{<:RealU} (nx) pixel centers
  • yc::Vector{<:RealU} (ny) pixel centers
  • mask::Matrix{Bool} (nx, ny) which pixels to reconstruct

option

  • ia_skip::Int default 1
  • T::Type{<:Number} usually same as eltype(sino)

out

  • image::Matrix (nx, ny)
source
Sinograms.fbp_back_par_xyMethod
fbp_back_par_xy(sino, sinϕ, cosϕ,
    wb, x, y ; T)

Pixel-driven back-projection for a single (x,y) location for sinogram sino from a parallel-beam geometry. It assumes the angles are equally spaced over [0,π).

in

  • sino::Matrix{<:Number} (nb, na) usually ramp-filtered
  • sinϕ::Vector{<:Real} (na)
  • cosϕ::Vector{<:Real} (na)
  • wb::Real = (nb+1)/2 + offset where usually offset=0
  • x,y::Real pixel center location, normalized by ray spacing

option

  • T::Type{<:Number} typically same as eltype(sino)

out

  • Returns a scalar of type T.
source
Sinograms.fbp_filterMethod
Hk = fbp_filter(rg::RayGeom ;
    npad=0, ds::RealU = rg.d, decon1::Bool=true, window=Window())

Compute frequency response of ramp-like filter used for FBP image reconstruction. Supports parallel-beam and fan-beam tomographic geometries in 2D and 3D. This code samples the band-limited ramp to avoid the aliasing that would be caused by sampling the ramp directly in the frequency domain.

in

  • rg::RayGeom

option

  • npad::Int # of padded samples. (default: next power of 2)
  • ds::Td detector sample spacing (default from st)
  • decon1::Bool deconvolve effect of linear interpolator? (default: true)
  • window::Window apodizer; default: Window()

out

  • Hk::Vector apodized ramp filter frequency response
source
Sinograms.fbp_rampFunction
h, n = fbp_ramp(rg::SinoGeom, N::Int)

'ramp-like' filters for parallel-beam and fan-beam FBP reconstruction. This sampled band-limited approach avoids the aliasing that would be caused by sampling the ramp directly in the frequency domain.

in

  • rg::SinoGeom
  • N::Int : # of samples (must be even)

out

  • h::Vector{<:RealU} : samples of band-limited ramp filter
  • n::UnitRange{Int64} : -(N÷2):(N÷2-1)
source
Sinograms.fbp_sino_filterMethod
sino = fbp_sino_filter(sino::Array, filter::Vector ; extra=0)

Apply ramp-like filters to sinogram(s) for 2D FBP image reconstruction. Supports both parallel-beam and fan-beam tomographic geometries in 2D and 3D.

in

  • sino::AbstractArray{<:Number} [nb (L)] sinograms
  • filter::AbstractVector (npad ≥ nb) apodized ramp filter frequency response

option

  • extra::Int # of extra sinogram radial samples to keep (default: 0)
  • npad::Int # of padded samples. (default: next power of 2)

out

  • sino::AbstractArray sinogram with filtered rows
source
Sinograms.fbp_windowMethod
fbp_window(w::Window, N::Int ; T = Float32)

Create an apodizing window of length N and fftshift it.

julia> fbp_window(Window(Hamming()), 4)
4-element Vector{Float32}:
 1.0
 0.54
 0.0
 0.54
source
Sinograms.fdkMethod
image = fdk(plan, proj)

Reconstruct 3D image from cone-beam computed tomography (CBCT) data collected with a circular source trajectory via FDK method.

in

  • plan::FDKplan
  • proj (ns,nt,na) projection views

out

  • image (nx,ny,nz) reconstructed image

References: Feldkamp, Davis, Kress, JOSA-A, 1(6):612-9, June 1984.

source
Sinograms.fdk_weight_cylFunction
fdk_weight_cyl

FDK projection weighting providing "exact" CBCT reconstruction for cylindrical-like objects that satisfy $f(x,y,z) = f(x,y,0) ∀z$. The output is a (ns,nt) matrix.

source
Sinograms.fft_filterFunction
fft_filter(data::Array, filter::Vector [, dim])

Apply filter to selected dimensions of array data using FFT.

in

  • data::AbstractArray{<:Number} (n, (L))
  • filter::AbstractVector (n) apodized ramp filter frequency response

option

  • dim non-singleton dimensions of filter (typically 1)

out

  • out::AbstractArray data filtered along dimension dim

If the input data is real, so will be the output; this assumes filter frequency response has appropriate symmetry.

source
Sinograms.footprint_sizeFunction
footprint_size(st::Union{SinoGeom,CtGeom}, ig::ImageGeom)

Unitless maximum footprint size (in detector pixels).

source
Sinograms.integrate1DMethod
integrate1D(p_v::Vector, pixelSize)

Calculates the integral image set, for a given column vector and a pixel size

source
Sinograms.map2xMethod
map2x(x1,y1,x2,y2)

Maps detector or pixel boundaries onto x-axis, for tube and detector rotation angles and detector/pixel boundaries

source
Sinograms.map2yMethod
map2y(x1,y1,x2,y2)

Maps detector or pixel boundaries onto y-axis, for tube and detector rotation angles and detector/pixel boundaries

source
Sinograms.oversampleMethod
oversample(rg, over::Int)

Over-sample sinogram geometry in "radial" dimension. For Mojette sampling, it means that d = dx/over.

source
Sinograms.parker_weightFunction
parker_weight(rg::SinoGeom, T = Float32)

Compute Parker weighting for non-360° orbits. See https://doi.org/10.1118/1.595078. Returns Matrix{T} of size:

  • (1,1) for SinoPar with typical 180 or 360 orbit
  • (1,na) for SinoPar with atypical orbit
  • (1,1) for SinoFan with typical 360 orbit
  • (ns,na) for SinoFan with typical 360 orbit
source
Sinograms.parker_weightMethod
parker_weight(rg::CtFan; T::Type{<:AbstractFloat} = Float32, kwargs...)

For 3D case, return Array{T,3} where size is

  • (1,1,1) typical fan case with 360° orbit
  • (ns,1,na) atypical fan case including short scan
source
Sinograms.plan_fbpMethod
plan = plan_fbp(rg, ig; window=Window(), ...)

Plan FDK 3D CBCT image reconstruction, with either flat or arc detector.

To use this method, you first call it with the CT geometry and image geometry. The routine returns the initialized plan. Thereafter, to to perform FDK reconstruction, call fbp with the plan (perhaps numerous times for the same geometry).

in

  • rg::CtGeom
  • ig::ImageGeom only reconstruct pixels within ig.mask.

option

  • window::Window e.g., Window(Hamming(), 0.8); default Window()
  • npad::Int # of radial bins after padding; default nextpow(2, rg.ns + 1)
  • decon1::Bool deconvolve interpolator effect? (default true)

out

  • plan::FDKplan initialized plan
source
Sinograms.plan_fbpMethod
plan = plan_fbp(rg, ig; how=:normal, window=Window())

Plan FBP 2D tomographic image reconstruction for parallel-beam & fan-beam cases, with either flat or arc detector for fan-beam case.

To use this method, you first call it with the sinogram geometry and image geometry. The routine returns the initialized plan. Thereafter, to to perform FBP reconstruction, call fbp with the plan (perhaps numerous times for the same geometry).

in

  • rg::SinoGeom
  • ig::ImageGeom only reconstruct pixels within ig.mask.

option

  • how::Symbol how to reconstruct
    • :normal default
    • :mojette use mojette rebinning and Gtomo2_table
  • window::Window e.g., Window(Hamming(), 0.8); default Window()
  • npad::Int # of radial bins after padding; default nextpow(2, rg.nb + 1)
  • decon1::Bool deconvolve interpolator effect? (default true)
  • T::Type{<:Number} type of sino elements (default Float32)

out

  • plan::FBPplan initialized plan
source
Sinograms.project_bddMethod
project_bdd(phantom::AbstractMatrix{<:T}, geo)

Generates a sinogram using the forward projection algorithm, for a phantom image and a tuple of geometry definitions

source
Sinograms.ramp_arcMethod
(h, n) = ramp_arc(N::Int, ds::RealU, dsd::RealU)

Ramp filter samples for arc fan geometry, for n = -(N÷2):(N÷2-1).

  • N must be even.
source
Sinograms.ramp_flatMethod
(h, n) = ramp_flat(N::Int, ds::RealU)

Ramp filter samples for flat fan geometry, for n = -(N÷2):(N÷2-1).

  • N must be even.
source
Sinograms.raysMethod
i = rays(rg::CtGeom)

Return parallel-beam coordinates of all rays for this CT geometry. Return type of i is a ProductIterator that makes tuples of the form (u, v, ϕ, θ). To make projections call p = [fun(c...) for c in i] where fun is radon(...).

source
Sinograms.raysMethod
i = rays(rg::SinoGeom)

Radial r and angular ϕ coordinates (in radians) of all sinogram elements for the given geometry. Return type of i is a ProductIterator that makes tuples of the form (r, ϕ). To make projections call p = [fun(c...) for c in i] where fun is radon(...).

source
Sinograms.to_radiansMethod
to_radians(angles::AbstractArray{<:AbstractFloat})

When Unitful package not loaded, assume angles are in degrees and convert to radians.

source
Sinograms.zwart_powellMethod
output = zwart_powell(r, ϕ)

Analytic 2D Radon transform value of Zwart-Powell box spline, for radial distance r (normalized by pixel size) and angle ϕ (in radians).

source