Methods list
ImagePhantoms.AbstractShape
ImagePhantoms.Cone
ImagePhantoms.Cuboid
ImagePhantoms.Cylinder
ImagePhantoms.Dirac
ImagePhantoms.Ellipse
ImagePhantoms.EllipsePhantomVersion
ImagePhantoms.Ellipsoid
ImagePhantoms.Gauss2
ImagePhantoms.Gauss3
ImagePhantoms.Object
ImagePhantoms.Object
ImagePhantoms.Object
ImagePhantoms.Object
ImagePhantoms.Object
ImagePhantoms.Object2d
ImagePhantoms.Object3d
ImagePhantoms.Rect
ImagePhantoms.SheppLogan
ImagePhantoms.SheppLoganToft
Base.:*
Base.show
ImagePhantoms.Rxyz_inv
ImagePhantoms.Rxyz_mul
ImagePhantoms._interval
ImagePhantoms.circle
ImagePhantoms.cone
ImagePhantoms.coords
ImagePhantoms.coords
ImagePhantoms.cube
ImagePhantoms.cube_bounds
ImagePhantoms.cuboid
ImagePhantoms.cylinder
ImagePhantoms.dirac2
ImagePhantoms.dirac3
ImagePhantoms.disk_phantom_params
ImagePhantoms.ellipse
ImagePhantoms.ellipse
ImagePhantoms.ellipse_parameters
ImagePhantoms.ellipse_parameters
ImagePhantoms.ellipse_parameters_shepplogan
ImagePhantoms.ellipse_parameters_uscale
ImagePhantoms.ellipsoid
ImagePhantoms.ellipsoid
ImagePhantoms.ellipsoid_parameters
ImagePhantoms.ellipsoid_parameters_shepplogan
ImagePhantoms.ellipsoid_parameters_uscale
ImagePhantoms.focus_chart
ImagePhantoms.fwhm2spread
ImagePhantoms.gauss2
ImagePhantoms.gauss3
ImagePhantoms.jinc
ImagePhantoms.mri_smap_basis
ImagePhantoms.mri_smap_fit
ImagePhantoms.mri_spectra
ImagePhantoms.phantom
ImagePhantoms.phantom
ImagePhantoms.phantom
ImagePhantoms.phantom
ImagePhantoms.phantom
ImagePhantoms.phantom
ImagePhantoms.phantom
ImagePhantoms.phantom
ImagePhantoms.phantom
ImagePhantoms.phantom1
ImagePhantoms.phantom1
ImagePhantoms.phantom1
ImagePhantoms.phantom1
ImagePhantoms.phantom1
ImagePhantoms.phantom1
ImagePhantoms.phantom1
ImagePhantoms.phantom1
ImagePhantoms.phantom1
ImagePhantoms.radon
ImagePhantoms.radon
ImagePhantoms.radon
ImagePhantoms.radon
ImagePhantoms.radon
ImagePhantoms.radon
ImagePhantoms.radon
ImagePhantoms.radon
ImagePhantoms.radon
ImagePhantoms.radon_tri
ImagePhantoms.radon_type
ImagePhantoms.rect
ImagePhantoms.rotate
ImagePhantoms.rotate
ImagePhantoms.scale
ImagePhantoms.shepp_logan
ImagePhantoms.shepp_logan_values
ImagePhantoms.spectrum
ImagePhantoms.spectrum
ImagePhantoms.spectrum
ImagePhantoms.spectrum
ImagePhantoms.spectrum
ImagePhantoms.spectrum
ImagePhantoms.spectrum
ImagePhantoms.spectrum
ImagePhantoms.spectrum
ImagePhantoms.spectrum1
ImagePhantoms.spectrum1
ImagePhantoms.spectrum1
ImagePhantoms.spectrum1
ImagePhantoms.spectrum1
ImagePhantoms.spectrum1
ImagePhantoms.spectrum1
ImagePhantoms.spectrum1
ImagePhantoms.spectrum1
ImagePhantoms.spectrum1
ImagePhantoms.sphere
ImagePhantoms.sphere_transform
ImagePhantoms.square
ImagePhantoms.translate
ImagePhantoms.trapezoid
ImagePhantoms.triangle
Methods usage
ImagePhantoms.AbstractShape
— TypeAbstractShape{D}
Generic shape type for ImagePhantoms
. The dimension D
is likely 2
or 3
, e.g., 2
for an Ellipse
and 3
for an Ellipsoid
.
ImagePhantoms.Cone
— TypeCone <: AbstractShape{3}
ImagePhantoms.Cuboid
— TypeCuboid <: AbstractShape{3}
ImagePhantoms.Cylinder
— TypeCylinder <: AbstractShape{3}
ImagePhantoms.Dirac
— TypeDirac{D} <: AbstractShape{D}
ImagePhantoms.Ellipse
— TypeEllipse <: AbstractShape{2}
ImagePhantoms.EllipsePhantomVersion
— TypeEllipsePhantomVersion
Parent type for different versions of ellipse phantoms:
SheppLogan
original CT version from Shepp&Logan paperSheppLoganEmis
higher contrast version suitable for emission tomographySheppLoganBrainWeb
integer index version based on brainwebSheppLoganToft
higher contrast version from Toft, 1996SouthPark
for fun
ImagePhantoms.Ellipsoid
— TypeEllipsoid <: AbstractShape{3}
ImagePhantoms.Gauss2
— TypeGauss2 <: AbstractShape{2}
ImagePhantoms.Gauss3
— TypeGauss3 <: AbstractShape{3}
ImagePhantoms.Object
— TypeObject(shape ; cx, cy, wx=1, wy=wx, ϕ=0, value=1)
2D object constructor from values (without tuples).
ImagePhantoms.Object
— TypeObject{S, D, V, C, A, Da}(center, width, angle, value) <: AbstractObject
where {S <: AbstractShape, V <: Number, C,A <: RealU}
General container for 2D and 3D objects for defining image phantoms.
center::NTuple{D,C}
coordinates of "center" of this objectwidth::NTuple{D,C}
"width" along each axis (e.g., FWHM for Gauss, radii for Ellipse)angle::NTuple{D-1,A}
angle of x' axis relative to x axis, in radians (or with units)value::V
"intensity" value for this object
Example
julia> Object(Ellipse(), (0,0), (1,2), 0.0, 1//2)
Object2d{Ellipse, Rational{Int64}, Int64, Float64, 1, Float64} (S, D, V, ...)
center::NTuple{2,Int64} (0, 0)
width::NTuple{2,Int64} (1, 2)
angle::Tuple{Float64} (0.0,)
value::Rational{Int64} 1//2
ImagePhantoms.Object
— TypeObject(shape ; cx, cy, cz, wx=1, wy=wx, wz=wx, ϕ=0, θ=0, ψ=0, value=1)
3D object constructor from values (without tuples).
ImagePhantoms.Object
— MethodObject(shape, center=(0,…), width=(1,…), angle=(0,…), value=1)
Object(shape ; center, width=(1,…), angle=(0,…), value=1)
General outer Object
constructor from tuples, as either positional arguments or named keyword arguments.
ImagePhantoms.Object
— MethodObject(ob::Object ; center, width, angle, value)
Make a copy of Object ob
, optionally modifying some values.
ImagePhantoms.Object2d
— TypeObject2d{S,V,C} = Object{S,2,V,C} where {S <: AbstractShape, V,C <: Number}
For 2D objects
ImagePhantoms.Object3d
— TypeObject3d{S,V,C} = Object{S,3,V,C} where {S <: AbstractShape, V,C <: Number}
For 3D objects
ImagePhantoms.Rect
— TypeRect <: AbstractShape{2}
ImagePhantoms.SheppLogan
— TypeSheppLogan
Original version from : Larry A Shepp, Benjamin F Logan, "The Fourier reconstruction of a head section," IEEE Transactions on Nuclear Science, 21(3):21-42, June 1974. doi
Also in Kak and Slaney 1988 text, p. 255. doi
ImagePhantoms.SheppLoganToft
— TypeSheppLoganToft
Toft, Peter Aundal & Sørensen, John Aasted "The Radon transform-theory and implementation," Technical University of Denmark (DTU), 1996. Page 201. https://files.openpdfs.org/1ra51GP6gJO.pdf
Base.:*
— Method(*)(ob::Object, x::Number)
(*)(x::Number, ob::Object)
Scale object value
by x
.
Base.show
— Methodshow(io::IO, ::MIME"text/plain", ob::Object)
ImagePhantoms.Rxyz_inv
— MethodRxyz_inv(x::RealU, y::RealU, z::RealU, ϕ::RealU, θ::RealU, ψ::RealU)
Rxyz_inv(x::RealU, y::RealU, sinϕ, sinθ, sinψ, cosϕ, cosθ, cosψ)
Multiply Rz(-ϕ) * Ry(-θ) * Rx(-ψ) * [x, y, z]
for 3D (inverse) rotation.
ImagePhantoms.Rxyz_mul
— MethodRxyz_mul(x::RealU, y::RealU, z::RealU, ϕ::RealU, θ::RealU, ψ::RealU)
Rxyz_mul(x::RealU, y::RealU, sinϕ, sinθ, sinψ, cosϕ, cosθ, cosψ)
Multiply Rx(ψ) * Ry(θ) * Rz(ϕ) * [x, y, z]
for 3D rotation, where ψ, θ, ϕ
denote rotation about x, y, z
axes with right-hand rule.
This is the reverse order of z-y′-x″ in https://en.wikipedia.org/wiki/Rotationformalismsinthreedimensions#Eulerangles(z-y%E2%80%B2-x%E2%80%B3intrinsic)%E2%86%92rotationmatrix.
ImagePhantoms._interval
— Method(lower, upper) = _interval(a, b)
Determine the interval [lower, upper]
corresponding to the set {x : b ≤ a x}
where a
is unitless but b
and x
have same units.
ImagePhantoms.circle
— Functioncircle(x, y, r, v=1) (circle of radius `r` centered at `(x,y)`)
circle((x,y), r=1, v=1) ditto
circle(r, v=1) centered at origin
Construct circles as special cases of Ellipse
.
ImagePhantoms.cone
— Methodcone(cx, cy, cz, wx, wy, wz, Φ, Θ, value::Number)
cone(center::NTuple{3,RealU}, width::NTuple{3,RealU}, angle::NTuple{3,RealU}, v)
Construct Object{Cone}
from parameters; here width
is the base radii for wx
and wy
and wz
is the height.
ImagePhantoms.coords
— Methodcoords(object::Object2d, x::RealU, y::RealU)
Put coordinates (x,y)
in canonical axes associated with object
.
ImagePhantoms.coords
— Methodcoords(object::Object3d, x::RealU, y::RealU, y::RealU)
Put coordinates (x,y,z)
in canonical axes associated with object
.
ImagePhantoms.cube
— Functioncube(x, y, z, w, v=1) (cube of width `w` centered at `(x,y,z)`)
cube((x,y,z), w, v=1) ditto
cube(w, v=1) centered at origin
Construct cubes as special cases of Cuboid
.
ImagePhantoms.cube_bounds
— Method(ℓmin, ℓmax) = cube_bounds(p, e)
Bounds of ℓ corresponding to rect(x) for direction e
from point p
.
ImagePhantoms.cuboid
— Methodcuboid(cx, cy, cz, wx, wy, wz, Φ, Θ, value::Number)
cuboid(center::NTuple{3,RealU}, width::NTuple{3,RealU}, angle::NTuple{3,RealU}, v)
Construct Object{Cuboid}
from parameters; here width
is the full-width.
ImagePhantoms.cylinder
— Methodcylinder(cx, cy, cz, wx, wy, wz, Φ, Θ, value::Number)
cylinder(center::NTuple{3,RealU}, width::NTuple{3,RealU}, angle::NTuple{3,RealU}, v)
Construct Object{Cylinder}
from parameters; here width
is the radius in x,y and the height in z.
ImagePhantoms.dirac2
— Methoddirac2(...)
Construct Object2d{Dirac2}
from parameters; here width
is a scale factor: δ((x - c)/w) = w * δ(x - c)
.
ImagePhantoms.dirac3
— Methoddirac3(...)
Construct Object3d{Dirac3}
from parameters; here width
is a scale factor: δ((x - c)/w) = w * δ(x - c)
.
ImagePhantoms.disk_phantom_params
— Methodparams = disk_phantom_params( ; ...)
Generate ndisk
ellipse phantom parameter 6-tuples for a head-sized disk plus many disks within it, designed so that the disks have some minimum separation minsep
to avoid overlap and to simplify patch-based model fitting.
Options
fov::Real = 240
image field of view in mmrhead::Function = () -> 100
background radius for "head" [mm]muhead::Function = () -> 1000
"μ" (intensity) value for background head diskrdisk::Function = () -> 10 + 10 * rand()
random disk radii [10,20]mudisk::Function = () -> 100 + 200 * rand()
"μ" values for disks [100,300]ndisk::Function = () -> 10
# of random disksminsep::Real = 8
minimum disk separation in mmmaxtry::Int = 500
give up on adding more disks if this is reachedwarn::Bool = false
warn if maxtry reached?seed::Int = 0
if nonzero then use this seed
The function options can be replaced with rand() for other Distributions.
ImagePhantoms.ellipse
— Methodellipse(cx, cy, rx=1, ry=rx, ϕ=0, value::Number=1)
ellipse(center::NTuple{2,RealU}, radii::NTuple{2,RealU}=(1,1), ϕ::RealU=0, v=1)
Construct Object{Ellipse}
from parameters.
ImagePhantoms.ellipse
— Methodphantom = ellipse(Vector{Tuple{6 params}})
Return vector of Object{Ellipse}
, one for each element of input vector of tuples. Often the input comes from ellipse_parameters
.
ImagePhantoms.ellipse_parameters
— Functionellipse_parameters(case; fovs::NTuple{2}, u::NTuple{3}, disjoint::Bool)
Return vector of Tuples of ellipse parameters. By default the first four elements of each tuple are unitless "fractions of field of view", so elements 1,3 are scaled by xfov
and elements 2,4 are scaled by yfov
, where (xfov, yfov) = fovs
. The optional 3-tuple u
specifies scaling and/or units:
- elements 1-4 (center, radii) are scaled by
u[1]
(e.g., mm), - elements 5 (angle) is scaled by
u[2]
(e.g.,1
or°
), - elements 6 (value) is scaled by
u[3]
(e.g.,1/cm
) for an attenuation map.
If disjoint==true
then the middle ellipse positions are adjusted to avoid overlap.
ImagePhantoms.ellipse_parameters
— Methodparams = ellipse_parameters(SouthPark() ; fovs::NTuple{2,Number} = (100,100))
Ellipse parameters for "South Park" phantom.
ImagePhantoms.ellipse_parameters_shepplogan
— Methodellipse_parameters_shepplogan( ; disjoint::Bool)
10 × 6 Matrix{Float64}
of classic Shepp-Logan ellipse parameters. If disjoint==true
then the middle ellipse positions are adjusted to avoid overlap.
ImagePhantoms.ellipse_parameters_uscale
— Methodellipse_parameters_uscale(params, fovs, uc, ua, uv)
Return vector of Tuples after FOV and unit scaling.
ImagePhantoms.ellipsoid
— Methodellipsoid(cx, cy, cz, rx=1, ry=1, rz=1, Φ=0, Θ=0, value::Number = 1)
ellipsoid(center::NTuple{3,RealU}, radii::NTuple{3,RealU}, angle::NTuple{3,RealU}, v)
Construct Object{Ellipsoid}
from parameters.
ImagePhantoms.ellipsoid
— Methodoa = ellipsoid(Vector{Tuple{10 params}})
Return vector of Object{Ellipsoid}
, one for each element of input vector of tuples. Often the input comes from ellipsoid_parameters
.
ImagePhantoms.ellipsoid_parameters
— Functionellipsoid_parameters(case; fovs::NTuple{3}, u::NTuple{3})
Return vector of Tuples of ellipsoid parameters. By default the parameters are for a 3D Shepp-Logan ellipsoid phantom. By default the first 6 elements of each tuple are unitless "fractions of field of view", so elements 1,4 are scaled by xfov
and elements 2,5 are scaled by yfov
, and elements 3,6 are scaled by zfov
, where (xfov, yfov, zfov) = fovs
. The optional 3-tuple u
specifies scaling and/or units:
- elements 1-6 (center, radii) are scaled by
u[1]
(e.g., mm), - elements 7-8 (angles) are scaled by
u[2]
(e.g.,1
or°
), - element 9 (value) is scaled by
u[3]
(e.g.,1/cm
) for an attenuation map.
ImagePhantoms.ellipsoid_parameters_shepplogan
— Methodellipsoid_parameters_shepplogan( ; disjoint::Bool)
12 × 10 Matrix{Float64}
of 3D Shepp-Logan ellipsoid parameters (cx,cy,cz, rx,ry rz, Φ,Θ=0,Ψ=0, ρ)
. By default the first 6 columns are unitless "fractions of field of view".
ImagePhantoms.ellipsoid_parameters_uscale
— Methodellipsoid_parameters_uscale(params, fovs, uc, ua, uv)
Return vector of Tuples after FOV and unit scaling.
ImagePhantoms.focus_chart
— Methodfocus_chart( ; radius=1, nspoke=30, value=1)
Generate nspoke
Triangle phantom parameters for a focus chart. Returns Vector{Object2d{Triangle}
for passing to phantom
.
Options
radius::RealU = 1
radius of phantomnspoke::Int = 60
# of spokesvalue::Number = 1
alternate between 0 and this value
ImagePhantoms.fwhm2spread
— Methods = fwhm2spread(w)
Convert FWHM w
to equivalent Gaussian spread s
for $\exp(-π (x/s)^2)$. exp(-π (fwhm/2/s)^2) = 1/2
means fwhm/2/s) = sqrt(log(2)/π)
2s/fwhm = sqrt(π/log(2))
2s = fwhm * sqrt(π/log(2))
s = fwhm * sqrt(π / log(16))
ImagePhantoms.gauss2
— Methodgauss2(cx, cy, wx, wy=wx, ϕ=0, value::Number=1)
gauss2(center::NTuple{2,RealU}, width::NTuple{2,RealU}=(1,1), ϕ::RealU=0, v=1)
gauss2(w, v=1) (isotropic of width `w`)
Construct Object{Gauss2}
from parameters; here width
= FWHM (full-width at half-maximum).
In 1D, the formula is g(x) = exp(-π ((x - cx) / sx)^2)
where sx = fwhm2spread(w) = w * sqrt(π / log(16))
, which, for cx=0
, has 1D FT G(ν) = sx^2 exp(π (sx νx)^2)
.
ImagePhantoms.gauss3
— Methodgauss3(cx, cy, cz, wx, wy, wz, Φ=0, Θ=0, value::Number = 1)
gauss3(center::NTuple{3,RealU}, radii::NTuple{3,RealU}=(1,1,1), angle::NTuple{2,RealU}=(0,0), v=1)
gauss3(r, v=1) (isotropic of width `w`)
Construct Object{Gauss3}
from parameters; here width
= FWHM (full-width at half-maximum).
ImagePhantoms.jinc
— Methodjinc(x)
Return jinc(x) = J1(π*x)/(2x)
, where J1
is a Bessel function of the first kind.
The argument x
must be unitless.
SpecialFunctions.jinc(0) = 1
whereas the convention here is jinc(0) = π / 4
which corresponds to the area of a disk of unit diameter.
ImagePhantoms.mri_smap_basis
— Methodmrismapbasis(mask ; kmax, kt, ki)
Construct Fourier basis for representing MRI sensitivity maps in terms of separable complex exponential signals, products of of the form basis = (k,N) -> exp.(2im * π * nfun(N) * kfun(k,N))
. The default is nfun(N) = -(N÷2):(N÷2)-1
which is suitable for even N
only. The default is kfun(k,N) = k / 2N
, which is DCT-II like frequencies, leading to better boundary behavior than the DFT frequencies k/N
.
Input
mask::AbstractArray{Bool,D}
binary support mask for region to reconstruct
Option
kmax::Int = 5
default frequency index -kmax:kmax in all dimensionskfun::Function = (k,N) -> k / (2N)
# DCT-II frequencydeltas::NTuple{D,<:Number} = ones(D)
pixel sizes
(For additional options kmaxs
, kt
, ki
, T
, see code.)
Output
(; B, ν)
whereB
is basis matrix of sizecount(mask) × nk
where typicallynk = (2*kmax+1)^D
andν
isnk
frequency tuples; each tuple has formν = kfun.(Tuple(k), size(mask)) ./ deltas
.
ImagePhantoms.mri_smap_fit
— Methodmri_smap_fit(smaps, embed ; mask, kwargs...)
Fit MRI sensitivity maps smaps
using mri_smap_basis(mask ; kwargs...)
. Caller provides ImageGeoms.embed
or equivalent.
Return named tuple (B, ν, coefs, nrmse, smaps)
:
(B, ν)
frommri_smap_basis
coefs::Vector
:[ncoil]
each of lengthnk
nrmse::Real
:smaps
vssmaps_fit
smaps::Vector{Array{D}}
:smaps_fit
ImagePhantoms.mri_spectra
— Methodmri_spectra(fx, fy, oa::Array{<:Object2d}, fit::NamedTuple)
Version of spectrum
suitable for parallel MRI with sensitivity maps that were fit previously using mri_smap_fit
. Returns a Vector of ncoil
kspace data Vectors of dimension length(fx)
ImagePhantoms.phantom
— Methodphantom(x::Vector, y::Vector, z::Vector, oa::Array{<:Object3d})
Return a 3D digital image of the phantom sampled at grid of (x,y,z)
locations. Returned 3D image size is length(x) × length(y) × length(z)
.
ImagePhantoms.phantom
— Methodphantom(x::Vector, y::Vector, oa::Array{<:Object2d}, oversample::Int; T)
Return a 2D digital image of the phantom sampled at grid of (x,y)
locations, with over-sampling factor oversample
and element type T
.
ImagePhantoms.phantom
— Methodphantom(x::Vector, y::Vector, oa::Array{<:Object2d})
Return a 2D digital image of the phantom sampled at grid of (x,y)
locations. Returned 2D image size is length(x) × length(y)
.
ImagePhantoms.phantom
— Methodphantom(itr, oa::Array{<:Object})
Return phantom values sampled at locations returned by generator (or iterator) itr
. Returned array size matches size(itr)
.
ImagePhantoms.phantom
— Methodphantom(oa::Array{<:Object2d})::Function
Return function of (x,y)
that user can sample at any locations to make a 2D phantom image for an Array
of 2D objects.
ImagePhantoms.phantom
— Methodphantom(oa::Array{<:Object3d})::Function
Return function of (x,y,z)
that user can sample at any locations to make a 3D phantom image for an Array
of 3D objects.
ImagePhantoms.phantom
— Methodphantom(ob::Object2d)::Function
Return function of (x,y)
that user can sample at any locations to make a 2D phantom image for a single 2D object.
ImagePhantoms.phantom
— Methodphantom(ob::Object3d)::Function
Return function of (x,y,z)
that user can sample at any locations to make a 3D phantom image for a single 3D object.
ImagePhantoms.phantom
— Methodimage = phantom(x, y, z, oa::Array{<:Object3d}, oversample::Int; T)
Return a 3D digital image of the phantom sampled at grid of (x,y,z)
locations, with over-sampling factor oversample
and element type T
.
ImagePhantoms.phantom1
— Methodphantom1(ob::Object2d{Ellipse}, (x,y))
Evaluate unit circle at (x,y)
, for unitless coordinates.
ImagePhantoms.phantom1
— Methodphantom1(ob::Object2d{Gauss2}, (x,y))
Evaluate unit gauss2 at (x,y)
, for unitless coordinates.
ImagePhantoms.phantom1
— Methodphantom1(ob::Object2d{Rect}, (x,y))
Evaluate unit square at (x,y)
, for unitless coordinates.
ImagePhantoms.phantom1
— Methodphantom1(ob::Object2d{Triangle}, (x,y))
Evaluate unit triangle at (x,y)
, for unitless coordinates.
ImagePhantoms.phantom1
— Methodphantom1(ob::Object3d{Cone}, (x,y,z))
Evaluate unit cone at (x,y,z)
, for unitless coordinates.
ImagePhantoms.phantom1
— Methodphantom1(ob::Object3d{Cuboid}, (x,y,z))
Evaluate unit cube at (x,y,z)
, for unitless coordinates.
ImagePhantoms.phantom1
— Methodphantom1(ob::Object3d{Cylinder}, (x,y,z))
Evaluate unit cylinder at (x,y,z)
, for unitless coordinates.
ImagePhantoms.phantom1
— Methodphantom1(ob::Object3d{Ellipsoid}, (x,y,z))
Evaluate unit sphere at (x,y,z)
, for unitless coordinates.
ImagePhantoms.phantom1
— Methodphantom1(ob::Object3d{Gauss3}, (x,y,z))
Evaluate Gauss3 (x,y,z)
, for unitless coordinates.
ImagePhantoms.radon
— Methodradon(u:Vector, v:Vector, ϕ:Vector, θ:Vector, oa::Array{<:Object3d})
Return parallel-beam projections sampled at grid of (u,v,ϕ,θ)
locations. Returned array size is length(u) × length(v) × length(ϕ) × length(θ)
.
ImagePhantoms.radon
— Methodradon(r::Vector, ϕ::Vector, oa::Array{<:Object2d})
Return parallel-beam 2D sinogram sampled at grid of (r,ϕ)
locations. Returned array size is length(r) × length(ϕ)
.
ImagePhantoms.radon
— Methodradon(u:Vector, v:Vector, ϕ:RealU, θ:RealU, oa::Array{<:Object3d})
Return parallel-beam projection view sampled at grid of (u,v)
locations for a given (ϕ,θ)
pair. Returned array size is length(u) × length(v)
.
ImagePhantoms.radon
— Methodradon(r::Vector, ϕ::RealU, oa::Array{<:Object2d})
Return parallel-beam projection sampled at grid of r
locations for one ϕ
value. Returned array size is length(r)
.
ImagePhantoms.radon
— Methodradon(itr, oa::Array{<:Object})
Return parallel-beam projections sampled at locations returned by generator (or iterator) itr
. Returned array size matches size(itr)
.
ImagePhantoms.radon
— Methodradon(oa::Array{<:Object2d})::Function
Return function of (r,ϕ)
that user can sample at any sinogram coordinates to make a phantom 2D sinogram for an array of 3D objects.
ImagePhantoms.radon
— Methodradon(oa::Array{<:Object3d})::Function
Return function of (u,v,ϕ,θ)
that user can sample at any projection coordinates to make projection views for an array of 3D objects.
ImagePhantoms.radon
— Methodradon(ob::Object2d)::Function
Return function of (r,ϕ)
that user can sample at any projection coordinates to make sinogram views of a 2D object.
The coordinate system used here is such that ϕ=0
corresponds to line integrals along the $y$ axis for an object $f(x,y)$. Then as ϕ
increases, the line integrals rotate counter-clockwise.
ImagePhantoms.radon
— Methodradon(ob::Object3d)::Function
Return function of (u,v,ϕ,θ)
that user can sample at any projection coordinates to make projection views of a 3D object.
The coordinate system used here is such that ϕ=0
corresponds to line integrals along the $y$ axis for an object $f(x,y,z)$. Then as ϕ
increases, the line integrals rotate counter-clockwise.
ImagePhantoms.radon_tri
— Methodradon_tri(r, sinϕ, cosϕ)
For a line integral at radial position r
and angle ϕ
, the locus of points along the line is {(r cos(ϕ), r sin(ϕ)) + τ (-sin(ϕ), cos(ϕ)) : τ ∈ ℝ}
. This function treats the equilateral triangle with base [-1/2,1/2], pointing upwards, as the intersection of three half planes:
- H0 =
{(x,y) : y ≥ 0}
- H1 =
y ≤ √3 (1/2 - x)
- H2 =
y ≤ √3 (1/2 + x)
.
Find the τ
values where the line locus lies in each half planes, then take the length of the intersection of those three intervals.
For example, for H1 we have r sin(ϕ) + τ cos(ϕ) ≤ √3 (1/2 - (r cos(ϕ) - τ sin(ϕ)))
or equivalently r (sin(ϕ) + √3 cos(ϕ)) - √3/2 ≤ τ (√3 sin(ϕ) - cos(ϕ))
which is a set the form b1 ≤ a1 τ
, corresponding to some interval (l1,u1)
. Similarly for H0 and H2.
This approach might not be the most efficient, but it is simple.
See Peter Aundal Toft, "The Radon transform - theory and implementation", 1996 https://orbit.dtu.dk/en/publications/the-radon-transform-theory-and-implementation for a different approach to finding the Radon transform of a triangle.
ImagePhantoms.radon_type
— Methodradon_type(::Object)
Determine the element type of the Radon transform of an object (including units if applicable). Ensures that its precision is at least Float32
.
ImagePhantoms.rect
— Methodrect(cx, cy, wx=1, wy=wx, ϕ=0, value::Number=1)
rect(center::NTuple{2,RealU}, width::NTuple{2,RealU}=(1,1), ϕ::RealU=0, v=1)
Construct Object{Rect}
from parameters; here width
is the full-width.
ImagePhantoms.rotate
— Methodrotate(ob::Object2d, θ::RealU)
Rotate a 2D object.
ImagePhantoms.rotate
— Methodrotate(ob::Object3d, (α,β))
rotate(ob::Object3d, α, β=0)
Rotate a 3D object.
ImagePhantoms.scale
— Methodscale(ob::Object, factor::RealU)
scale(ob::Object, factor::NTuple{D,RealU})
Scale the width(s) by factor
.
ImagePhantoms.shepp_logan
— Methodimage = shepp_logan(M, [N,] [case] ; options...)
Convenience method for generating M×N
samples of Shepp-Logan phantoms.
In
M::Int
: horizontal sizeN::Int
: vertical size, defaults toM
case::EllipsePhantomVersion = SheppLogan()
Options
oversample::Int = 3
(usually)yflip::Bool = true
(reversey
samples for convenience.)T::Type{<:Number}
defaultFloat32
(exceptInt
forBrainWeb
version)kwargs...
remaining options passed toellipse_parameters
for parameters.
Out
image
:M × N
matrix
The default here is 3× over-sampling along both axes (9 samples per pixel), except for the SheppLoganBrainWeb
phantom that consists of integer indices.
ImagePhantoms.shepp_logan_values
— Methodvalues = shepp_logan_values(::EllipsePhantomVersion)
Return 10 Shepp-Logan ellipse amplitudes for various versions.
ImagePhantoms.spectrum
— Methodspectrum(fx::Vector, fy::Vector, fz::Vector, oa::Array{<:Object3d})
Return 3D k-space array sampled at grid of (fx,fy,fz)
locations. Returned 3D array size is length(fx) × length(fy) × length(fz)
.
ImagePhantoms.spectrum
— Methodspectrum(fx:Vector, fy:Vector, oa::Array{<:Object2d})
Return 2D k-space array sampled at grid of (fx,fy)
locations. Returned 2D array size is length(fx) × length(fy)
.
ImagePhantoms.spectrum
— Methodspectrum(itr, oa::Array{<:Object})
Return spectrum of object(s) sampled at k-space locations returned by generator (or iterator) itr
. Returned array size matches size(itr)
.
ImagePhantoms.spectrum
— Methodspectrum(oa::Array{<:Object2d})::Function
Return function kspace(fx,fy)
that user can sample at any spatial frequency locations to evaluate the spectrum (2D Fourier transform) for an Array
of 2D objects.
ImagePhantoms.spectrum
— Methodspectrum(oa::Array{<:Object3d})::Function
Return function kspace(fx,fy,fz)
that user can sample at any spatial frequency locations to evaluate the spectrum (3D Fourier transform) for an Array
of 3D objects.
ImagePhantoms.spectrum
— Methodspectrum(ob::Object2d, coefs::AbstractVector, f::)
Version of spectrum(ob)
suitable for parallel MRI with sensitivity maps that were fit previously using mri_smap_fit
for a single coil with fit coefficients coefs
and frequencies f
(array of tuples).
ImagePhantoms.spectrum
— Methodspectrum(ob::Object2d, fit::NamedTuple, coil::Int)
Version of spectrum(ob)
suitable for parallel MRI with sensitivity maps that were fit previously using mri_smap_fit
.
ImagePhantoms.spectrum
— Methodspectrum(ob::Object2d)::Function
Return function of (fx,fy)
that user can sample at any spatial frequency locations to evaluate the spectrum (2D Fourier transform) of a single 2D object. Units of spatial frequencies should be the reciprocal of the units defining the object.
ImagePhantoms.spectrum
— Methodspectrum(ob::Object3d)::Function
Return function of (fx,fy,fz)
that user can sample at any spatial frequency locations to evaluate the spectrum (3D Fourier transform) of a single 3D object. Units of spatial frequencies should be the reciprocal of the units defining the object.
ImagePhantoms.spectrum1
— Methodspectrum1(::Object2d{Dirac2}, (kx,ky))
Spectrum of Dirac2 at (kx,ky)
, for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1
— Methodspectrum1(::Object2d{Ellipse}, (kx,ky))
Spectrum of unit circle at (kx,ky)
, for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1
— Methodspectrum1(::Object2d{Gauss2}, (kx,ky))
Spectrum of unit gauss2 at (kx,ky)
, for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1
— Methodspectrum1(::Object2d{Rect}, (kx,ky))
Spectrum of unit square at (kx,ky)
, for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1
— Methodspectrum1(ob::Object2d{Triangle}, (kx,ky))
Spectrum of unit triangle at (kx,ky)
, for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1
— Methodspectrum1(::Object3d{Cuboid}, (kx,ky,kz))
Spectrum of unit cube at (kx,ky,kz)
, for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1
— Methodspectrum(::Object3d{Cylinder}, (kx,ky,kz))
Spectrum of unit cylinder at (kx,ky,kz)
, for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1
— Methodspectrum1(::Object3d{Dirac3}, (kx,ky,kz))
Spectrum of Dirac3 at (kx,ky,kz)
, for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1
— Methodspectrum1(::Object3d{Ellipsoid}, (kx,ky,kz))
Spectrum of unit sphere at (kx,ky,kz)
, for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1
— Methodspectrum1(::Object3d{Gauss3}, (kx,ky,kz))
Spectrum of unit sphere at (kx,ky,kz)
, for unitless spatial frequency coordinates.
ImagePhantoms.sphere
— Functionsphere(x, y, z, r,v=1) (sphere of radius `r` centered at `(x,y,z)`)
sphere((x,y,z), r, v=1) ditto
sphere(r, v=1) centered at origin
Construct spheres as special cases of Ellipsoid
.
ImagePhantoms.sphere_transform
— Methodsphere_transform(f::Real)
Fourier transform of unit-radius sphere. The argument f
is the radial coordinate in k-space and is unitless. See p253 of Bracewell 1978, The Fourier transform and its applications, or https://doi.org/10.1002/mrm.21292.
Formula: 4/3 π
for f ≈ 0
, otherwise (sin(2πf) - 2πf cos(2πf)) / (2 * π^2 * f^3)
.
ImagePhantoms.square
— Functionsquare(x, y, w,v=1) (square of width `w` centered at `(x,y)`)
square((x,y), w=1, v=1) ditto
square(w, v=1) centered at origin
Construct squares as special cases of Rect
.
ImagePhantoms.translate
— Methodtranslate(ob::Object, shift::NTuple{D,RealU})
translate(ob::Object{S,2}, xshift, yshift)
translate(ob::Object{S,3}, xshift, yshift, zshift)
Translate the center coordinates of an object by shift
ImagePhantoms.trapezoid
— Methodtrapezoid(t::Real, t1, t2, t3, t4)
Unit-height trapezoid with breakpoints t1
, t2
, t3
, t4
.
ImagePhantoms.triangle
— Methodtriangle(cx, cy, wx=1, wy=wx, ϕ=0, value::Number=1, param::Real=0.5)
triangle(center::NTuple{2,RealU}, width::NTuple{2,RealU}=(1,1), ϕ::RealU=0, v=1, param=0.5)
Construct Object{Triangle}
from parameters. In the typical case where param=0.5
and width[1] == width[2]
, this is an equilateral triangle with base width[1]
centered along the x axis.