Methods list
ImagePhantoms.AbstractShapeImagePhantoms.ConeImagePhantoms.CuboidImagePhantoms.CylinderImagePhantoms.DiracImagePhantoms.EllipseImagePhantoms.EllipsePhantomVersionImagePhantoms.EllipsoidImagePhantoms.Gauss2ImagePhantoms.Gauss3ImagePhantoms.ObjectImagePhantoms.ObjectImagePhantoms.ObjectImagePhantoms.ObjectImagePhantoms.ObjectImagePhantoms.Object2dImagePhantoms.Object3dImagePhantoms.RectImagePhantoms.SheppLoganImagePhantoms.SheppLoganToftBase.:*Base.showImagePhantoms.Rxyz_invImagePhantoms.Rxyz_mulImagePhantoms._intervalImagePhantoms.circleImagePhantoms.coneImagePhantoms.coordsImagePhantoms.coordsImagePhantoms.cubeImagePhantoms.cube_boundsImagePhantoms.cuboidImagePhantoms.cylinderImagePhantoms.dirac2ImagePhantoms.dirac3ImagePhantoms.disk_phantom_paramsImagePhantoms.ellipseImagePhantoms.ellipseImagePhantoms.ellipse_parametersImagePhantoms.ellipse_parametersImagePhantoms.ellipse_parameters_shepploganImagePhantoms.ellipse_parameters_uscaleImagePhantoms.ellipsoidImagePhantoms.ellipsoidImagePhantoms.ellipsoid_parametersImagePhantoms.ellipsoid_parameters_shepploganImagePhantoms.ellipsoid_parameters_uscaleImagePhantoms.focus_chartImagePhantoms.fwhm2spreadImagePhantoms.gauss2ImagePhantoms.gauss3ImagePhantoms.jincImagePhantoms.mri_smap_basisImagePhantoms.mri_smap_fitImagePhantoms.mri_spectraImagePhantoms.phantomImagePhantoms.phantomImagePhantoms.phantomImagePhantoms.phantomImagePhantoms.phantomImagePhantoms.phantomImagePhantoms.phantomImagePhantoms.phantomImagePhantoms.phantomImagePhantoms.phantom1ImagePhantoms.phantom1ImagePhantoms.phantom1ImagePhantoms.phantom1ImagePhantoms.phantom1ImagePhantoms.phantom1ImagePhantoms.phantom1ImagePhantoms.phantom1ImagePhantoms.phantom1ImagePhantoms.radonImagePhantoms.radonImagePhantoms.radonImagePhantoms.radonImagePhantoms.radonImagePhantoms.radonImagePhantoms.radonImagePhantoms.radonImagePhantoms.radonImagePhantoms.radon_triImagePhantoms.radon_typeImagePhantoms.rectImagePhantoms.rotateImagePhantoms.rotateImagePhantoms.scaleImagePhantoms.shepp_loganImagePhantoms.shepp_logan_valuesImagePhantoms.spectrumImagePhantoms.spectrumImagePhantoms.spectrumImagePhantoms.spectrumImagePhantoms.spectrumImagePhantoms.spectrumImagePhantoms.spectrumImagePhantoms.spectrumImagePhantoms.spectrumImagePhantoms.spectrum1ImagePhantoms.spectrum1ImagePhantoms.spectrum1ImagePhantoms.spectrum1ImagePhantoms.spectrum1ImagePhantoms.spectrum1ImagePhantoms.spectrum1ImagePhantoms.spectrum1ImagePhantoms.spectrum1ImagePhantoms.spectrum1ImagePhantoms.sphereImagePhantoms.sphere_transformImagePhantoms.squareImagePhantoms.translateImagePhantoms.trapezoidImagePhantoms.triangle
Methods usage
ImagePhantoms.AbstractShape — Type
AbstractShape{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 — Type
Cone <: AbstractShape{3}ImagePhantoms.Cuboid — Type
Cuboid <: AbstractShape{3}ImagePhantoms.Cylinder — Type
Cylinder <: AbstractShape{3}ImagePhantoms.Dirac — Type
Dirac{D} <: AbstractShape{D}ImagePhantoms.Ellipse — Type
Ellipse <: AbstractShape{2}ImagePhantoms.EllipsePhantomVersion — Type
EllipsePhantomVersionParent type for different versions of ellipse phantoms:
SheppLoganoriginal CT version from Shepp&Logan paperSheppLoganEmishigher contrast version suitable for emission tomographySheppLoganBrainWebinteger index version based on brainwebSheppLoganTofthigher contrast version from Toft, 1996SouthParkfor fun
ImagePhantoms.Ellipsoid — Type
Ellipsoid <: AbstractShape{3}ImagePhantoms.Gauss2 — Type
Gauss2 <: AbstractShape{2}ImagePhantoms.Gauss3 — Type
Gauss3 <: AbstractShape{3}ImagePhantoms.Object — Type
Object(shape ; cx, cy, cz, wx=1, wy=wx, wz=wx, ϕ=0, θ=0, ψ=0, value=1)3D object constructor from values (without tuples).
ImagePhantoms.Object — Type
Object(shape ; cx, cy, wx=1, wy=wx, ϕ=0, value=1)2D object constructor from values (without tuples).
ImagePhantoms.Object — Type
Object{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//2ImagePhantoms.Object — Method
Object(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 — Method
Object(ob::Object ; center, width, angle, value)Make a copy of Object ob, optionally modifying some values.
ImagePhantoms.Object2d — Type
Object2d{S,V,C} = Object{S,2,V,C} where {S <: AbstractShape, V,C <: Number}For 2D objects
ImagePhantoms.Object3d — Type
Object3d{S,V,C} = Object{S,3,V,C} where {S <: AbstractShape, V,C <: Number}For 3D objects
ImagePhantoms.Rect — Type
Rect <: AbstractShape{2}ImagePhantoms.SheppLogan — Type
SheppLoganOriginal 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 — Type
SheppLoganToftToft, 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
ImagePhantoms.Rxyz_inv — Method
Rxyz_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 — Method
Rxyz_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 — Function
circle(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 originConstruct circles as special cases of Ellipse.
ImagePhantoms.cone — Method
cone(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 — Method
coords(object::Object2d, x::RealU, y::RealU)Put coordinates (x,y) in canonical axes associated with object.
ImagePhantoms.coords — Method
coords(object::Object3d, x::RealU, y::RealU, y::RealU)Put coordinates (x,y,z) in canonical axes associated with object.
ImagePhantoms.cube — Function
cube(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 originConstruct 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 — Method
cuboid(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 — Method
cylinder(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 — Method
dirac2(...)Construct Object2d{Dirac2} from parameters; here width is a scale factor: δ((x - c)/w) = w * δ(x - c).
ImagePhantoms.dirac3 — Method
dirac3(...)Construct Object3d{Dirac3} from parameters; here width is a scale factor: δ((x - c)/w) = w * δ(x - c).
ImagePhantoms.disk_phantom_params — Method
params = 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 = 240image field of view in mmrhead::Function = () -> 100background 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 = 8minimum disk separation in mmmaxtry::Int = 500give up on adding more disks if this is reachedwarn::Bool = falsewarn if maxtry reached?seed::Int = 0if nonzero then use this seed
The function options can be replaced with rand() for other Distributions.
ImagePhantoms.ellipse — Method
ellipse(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 — Method
phantom = 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 — Function
ellipse_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.,1or°), - 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 — Method
params = ellipse_parameters(SouthPark() ; fovs::NTuple{2,Number} = (100,100))Ellipse parameters for "South Park" phantom.
ImagePhantoms.ellipse_parameters_shepplogan — Method
ellipse_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 — Method
ellipse_parameters_uscale(params, fovs, uc, ua, uv)Return vector of Tuples after FOV and unit scaling.
ImagePhantoms.ellipsoid — Method
ellipsoid(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 — Method
oa = 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 — Function
ellipsoid_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.,1or°), - element 9 (value) is scaled by
u[3](e.g.,1/cm) for an attenuation map.
ImagePhantoms.ellipsoid_parameters_shepplogan — Method
ellipsoid_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 — Method
ellipsoid_parameters_uscale(params, fovs, uc, ua, uv)Return vector of Tuples after FOV and unit scaling.
ImagePhantoms.focus_chart — Method
focus_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 = 1radius of phantomnspoke::Int = 60# of spokesvalue::Number = 1alternate between 0 and this value
ImagePhantoms.fwhm2spread — Method
s = 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 — Method
gauss2(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 — Method
gauss3(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 — Method
jinc(x)Return jinc(x) = J1(π*x)/(2x), where J1 is a Bessel function of the first kind.
The argument x must be unitless.
ImagePhantoms.mri_smap_basis — Method
mrismapbasis(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 = 5default 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, ν)whereBis basis matrix of sizecount(mask) × nkwhere typicallynk = (2*kmax+1)^Dandνisnkfrequency tuples; each tuple has formν = kfun.(Tuple(k), size(mask)) ./ deltas.
ImagePhantoms.mri_smap_fit — Method
mri_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_basiscoefs::Vector:[ncoil]each of lengthnknrmse::Real:smapsvssmaps_fitsmaps::Vector{Array{D}}:smaps_fit
ImagePhantoms.mri_spectra — Method
mri_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 — Method
phantom(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 — Method
phantom(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 — Method
phantom(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 — Method
phantom(itr, oa::Array{<:Object})Return phantom values sampled at locations returned by generator (or iterator) itr. Returned array size matches size(itr).
ImagePhantoms.phantom — Method
phantom(oa::Array{<:Object2d})::FunctionReturn 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 — Method
phantom(oa::Array{<:Object3d})::FunctionReturn 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 — Method
phantom(ob::Object2d)::FunctionReturn function of (x,y) that user can sample at any locations to make a 2D phantom image for a single 2D object.
ImagePhantoms.phantom — Method
phantom(ob::Object3d)::FunctionReturn 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 — Method
image = 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 — Method
phantom1(ob::Object2d{Ellipse}, (x,y))Evaluate unit circle at (x,y), for unitless coordinates.
ImagePhantoms.phantom1 — Method
phantom1(ob::Object2d{Gauss2}, (x,y))Evaluate unit gauss2 at (x,y), for unitless coordinates.
ImagePhantoms.phantom1 — Method
phantom1(ob::Object2d{Rect}, (x,y))Evaluate unit square at (x,y), for unitless coordinates.
ImagePhantoms.phantom1 — Method
phantom1(ob::Object2d{Triangle}, (x,y))Evaluate unit triangle at (x,y), for unitless coordinates.
ImagePhantoms.phantom1 — Method
phantom1(ob::Object3d{Cone}, (x,y,z))Evaluate unit cone at (x,y,z), for unitless coordinates.
ImagePhantoms.phantom1 — Method
phantom1(ob::Object3d{Cuboid}, (x,y,z))Evaluate unit cube at (x,y,z), for unitless coordinates.
ImagePhantoms.phantom1 — Method
phantom1(ob::Object3d{Cylinder}, (x,y,z))Evaluate unit cylinder at (x,y,z), for unitless coordinates.
ImagePhantoms.phantom1 — Method
phantom1(ob::Object3d{Ellipsoid}, (x,y,z))Evaluate unit sphere at (x,y,z), for unitless coordinates.
ImagePhantoms.phantom1 — Method
phantom1(ob::Object3d{Gauss3}, (x,y,z))Evaluate Gauss3 (x,y,z), for unitless coordinates.
ImagePhantoms.radon — Method
radon(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 — Method
radon(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 — Method
radon(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 — Method
radon(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 — Method
radon(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 — Method
radon(oa::Array{<:Object2d})::FunctionReturn 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 — Method
radon(oa::Array{<:Object3d})::FunctionReturn function of (u,v,ϕ,θ) that user can sample at any projection coordinates to make projection views for an array of 3D objects.
ImagePhantoms.radon — Method
radon(ob::Object2d)::FunctionReturn 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 — Method
radon(ob::Object3d)::FunctionReturn 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 — Method
radon_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 — Method
radon_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 — Method
rect(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 — Method
rotate(ob::Object2d, θ::RealU)Rotate a 2D object.
ImagePhantoms.rotate — Method
rotate(ob::Object3d, (α,β))
rotate(ob::Object3d, α, β=0)Rotate a 3D object.
ImagePhantoms.scale — Method
scale(ob::Object, factor::RealU)
scale(ob::Object, factor::NTuple{D,RealU})Scale the width(s) by factor.
ImagePhantoms.shepp_logan — Method
image = 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 toMcase::EllipsePhantomVersion = SheppLogan()
Options
oversample::Int = 3(usually)yflip::Bool = true(reverseysamples for convenience.)T::Type{<:Number}defaultFloat32(exceptIntforBrainWebversion)kwargs...remaining options passed toellipse_parametersfor parameters.
Out
image:M × Nmatrix
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 — Method
values = shepp_logan_values(::EllipsePhantomVersion)Return 10 Shepp-Logan ellipse amplitudes for various versions.
ImagePhantoms.spectrum — Method
spectrum(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 — Method
spectrum(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 — Method
spectrum(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 — Method
spectrum(oa::Array{<:Object2d})::FunctionReturn 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 — Method
spectrum(oa::Array{<:Object3d})::FunctionReturn 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 — Method
spectrum(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 — Method
spectrum(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 — Method
spectrum(ob::Object2d)::FunctionReturn 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 — Method
spectrum(ob::Object3d)::FunctionReturn 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 — Method
spectrum1(::Object2d{Dirac2}, (kx,ky))Spectrum of Dirac2 at (kx,ky), for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1 — Method
spectrum1(::Object2d{Ellipse}, (kx,ky))Spectrum of unit circle at (kx,ky), for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1 — Method
spectrum1(::Object2d{Gauss2}, (kx,ky))Spectrum of unit gauss2 at (kx,ky), for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1 — Method
spectrum1(::Object2d{Rect}, (kx,ky))Spectrum of unit square at (kx,ky), for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1 — Method
spectrum1(ob::Object2d{Triangle}, (kx,ky))Spectrum of unit triangle at (kx,ky), for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1 — Method
spectrum1(::Object3d{Cuboid}, (kx,ky,kz))Spectrum of unit cube at (kx,ky,kz), for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1 — Method
spectrum(::Object3d{Cylinder}, (kx,ky,kz))Spectrum of unit cylinder at (kx,ky,kz), for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1 — Method
spectrum1(::Object3d{Dirac3}, (kx,ky,kz))Spectrum of Dirac3 at (kx,ky,kz), for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1 — Method
spectrum1(::Object3d{Ellipsoid}, (kx,ky,kz))Spectrum of unit sphere at (kx,ky,kz), for unitless spatial frequency coordinates.
ImagePhantoms.spectrum1 — Method
spectrum1(::Object3d{Gauss3}, (kx,ky,kz))Spectrum of unit sphere at (kx,ky,kz), for unitless spatial frequency coordinates.
ImagePhantoms.sphere — Function
sphere(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 originConstruct spheres as special cases of Ellipsoid.
ImagePhantoms.sphere_transform — Method
sphere_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 — Function
square(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 originConstruct squares as special cases of Rect.
ImagePhantoms.translate — Method
translate(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 — Method
trapezoid(t::Real, t1, t2, t3, t4)Unit-height trapezoid with breakpoints t1, t2, t3, t4.
ImagePhantoms.triangle — Method
triangle(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.