API Reference
This page provides a comprehensive reference of the types and functions in MicroscopePSFs.jl.
Public API
These types and functions form the core public API of MicroscopePSFs.jl.
PSF Types
MicroscopePSFs.AbstractPSF
— TypeAbstractPSF
Abstract base type for all Point Spread Functions
2D PSF Models
MicroscopePSFs.GaussianPSF
— TypeGaussianPSF{T<:AbstractFloat} <: Abstract2DPSF{T}
Isotropic 2D Gaussian PSF.
Fields
σ
: Standard deviation in physical units (typically microns)
MicroscopePSFs.AiryPSF
— TypeAiryPSF{T<:AbstractFloat} <: Abstract2DPSF{T}
2D Airy pattern PSF using paraxial, scalar model.
Field amplitude for a circular aperture is given by: A(r) = ν/√(4π) * (2J₁(νr)/(νr)) where: ν = 2π*nₐ/λ J₁ is the Bessel function of first kind, order 1 r is the radial distance from optical axis
Fields
nₐ
: Numerical apertureλ
: Wavelength in micronsν
: Optical parameter = 2π*nₐ/λ
3D PSF Models
MicroscopePSFs.ScalarPSF
— TypeScalarPSF{T} <: Abstract3DPSF{T}
Scalar 3D PSF using explicit pupil function representation.
Fields
nₐ::T
: Numerical apertureλ::T
: Wavelength in micronsn::T
: Refractive indexpupil::PupilFunction{T}
: Complex pupil functionzernike_coeffs::Union{Nothing, ZernikeCoefficients{T}}
: Zernike coefficients used to create this PSF (if applicable)
Notes
Can be initialized with either a PupilFunction or ZernikeCoefficients using the ScalarPSF factory function.
MicroscopePSFs.VectorPSF
— TypeVectorPSF{T<:AbstractFloat} <: Abstract3DPSF{T}
Vector PSF model using explicit pupil function representation. Allows direct manipulation of pupil function for custom aberrations.
Fields
nₐ::T
: Numerical apertureλ::T
: Wavelength in micronsn_medium::T
: Sample medium refractive indexn_coverslip::T
: Cover slip refractive indexn_immersion::T
: Immersion medium refractive indexdipole::DipoleVector{T}
: Dipole orientationz_stage::T
: Distance the sample stage was moved away from the nominal focal plane at the coverslip (μm)vector_pupils::Vector{VectorPupilFunction{T}}
: Vector of pupil functions containing field components. For a single dipole, this contains one pupil. For a rotating dipole, it contains three pupils for x, y, and z orientations.base_pupil::Union{Nothing, PupilFunction{T}}
: Base pupil function representing system aberrationszernike_coeffs::Union{Nothing, ZernikeCoefficients{T}}
: Zernike coefficients used to create this PSF (if applicable)
PSF Acceleration
MicroscopePSFs.SplinePSF
— TypeSplinePSF{T<:AbstractFloat, IT<:AbstractInterpolation} <: AbstractPSF
A point spread function (PSF) represented as a B-spline interpolation.
Fields
spline
: The B-spline interpolation objectx_range
: Range of x-coordinates used for uniform grid interpolationy_range
: Range of y-coordinates used for uniform grid interpolationz_range
: Range of z-coordinates for 3D PSFs, ornothing
for 2D PSFsoriginal_grid
: Original grid data used to create the interpolationinterp_order
: Interpolation order used (0=constant, 1=linear, 3=cubic)
Notes
- Coordinates and ranges are in physical units (typically microns)
- PSF values are preserved from the original PSF that was sampled
- Full implementation is in spline_psf.jl
Core Interface Functions
MicroscopePSFs.amplitude
— Functionamplitude(p::PupilFunction, x::Real, y::Real, z::Real)
Calculate complex amplitude from pupil function integration.
Arguments
p::PupilFunction
: Pupil functionx::Real
: X position in μmy::Real
: Y position in μmz::Real
: Z position in μm
amplitude(psf::AbstractPSF, x::Real, y::Real)
Evaluate complex field amplitude at position (x,y) relative to PSF center.
Arguments
x, y
: Position in microns relative to PSF center
Returns
- Complex amplitude normalized such that |amplitude|² gives normalized intensity
Coordinates
Input coordinates (x,y) are in physical units (microns) relative to PSF center.
Examples
psf = AiryPSF(1.4, 0.532)
amp = amplitude(psf, 0.1, 0.2)
intensity = abs2(amp) # Convert to intensity
amplitude(psf::AbstractPSF, x::Real, y::Real, z::Real)
Evaluate complex field amplitude in 3D. z is axial distance from focus in microns.
amplitude(psf::ScalarPSF{T}, x::Real, y::Real, z::Real) where {T}
Calculate complex amplitude at a 3D position using scalar diffraction theory.
Arguments
psf
: ScalarPSF instancex, y, z
: Position in microns relative to PSF center
Returns
- Complex amplitude at the specified position
Notes
- Uses Fourier optics to propagate from pupil to image space
- Accounts for defocus and aberrations encoded in the pupil function
amplitude(psf::VectorPSF, x::Real, y::Real, z::Real)
Compute complex vector amplitude at given position.
Arguments
psf
: Vector PSF instancex, y
: Lateral position in microns relative to PSF centerz
: Axial position in microns representing depth above the coverslip
Returns
- Vector [Ex, Ey] of complex field amplitudes
Notes
- z coordinate represents the depth above the coverslip
- z_stage in the PSF indicates the distance the stage was moved away from the nominal focal plane
- Includes both UAF and SAF contributions automatically
- Not meaningful for rotating dipoles (throws an error if PSF has multiple pupils)
amplitude(psf::SplinePSF, x::Real, y::Real, z::Real)
Calculate the complex amplitude of the 3D PSF at position (x, y, z).
Arguments
psf
: SplinePSF instancex, y, z
: Coordinates in microns relative to PSF center
Returns
- Complex amplitude = sqrt(intensity) with zero phase
Notes
- Returns sqrt(intensity) as Complex to match the interface of other PSFs
- The SplinePSF does not model phase information
amplitude(psf::SplinePSF, x::Real, y::Real)
Calculate the complex amplitude of the PSF at position (x, y) with z=0.
Arguments
psf
: SplinePSF instancex, y
: Coordinates in microns relative to PSF center
Returns
- Complex amplitude = sqrt(intensity) with zero phase
MicroscopePSFs.integrate_pixels
— Functionintegrate_pixels(
psf::AbstractPSF,
pixel_edges_x::AbstractVector,
pixel_edges_y::AbstractVector,
emitter::AbstractEmitter;
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF intensity over camera pixels with optional support region optimization.
For each pixel in the specified region, numerically integrates the PSF intensity using the specified sampling density. Physical coordinates are relative to camera with (0,0) at top-left corner. Automatically handles z-coordinate if both PSF and emitter support it.
Arguments
psf::AbstractPSF
: Point spread function to integratepixel_edges_x::AbstractVector
: X-coordinate edges of pixels in micronspixel_edges_y::AbstractVector
: Y-coordinate edges of pixels in micronsemitter::AbstractEmitter
: Emitter with position in microns relative to camerasupport
: Region to calculate (default: Inf = full image)- If Real: radius in microns around emitter
- If Tuple: explicit (xmin, xmax, ymin, ymax) in microns
sampling::Integer=2
: Number of samples per pixel in each dimensionthreaded::Bool=true
: Whether to use multi-threading for integration Set to false when using with automatic differentiation frameworks
Returns
Matrix{T}
: Integrated intensities where T matches emitter.photons type- Array is indexed as [y,x] with [1,1] at top-left pixel
- Values represent actual photon counts in each pixel based on emitter's photon value
Examples
# Create pixel edges for a 20x20 camera with 100nm pixels
pixel_edges_x = pixel_edges_y = 0:0.1:2.0
emitter = Emitter2D(1.0, 1.0, 1000.0) # Emitter at (1μm, 1μm) with 1000 photons
psf = GaussianPSF(0.15) # σ = 150nm
# Calculate over full image
pixels = integrate_pixels(psf, pixel_edges_x, pixel_edges_y, emitter)
# Calculate only within a 0.5μm radius of the emitter
pixels_roi = integrate_pixels(psf, pixel_edges_x, pixel_edges_y, emitter, support=0.5)
See also: integrate_pixels_amplitude
, AbstractPSF
integrate_pixels(
psf::AbstractPSF,
camera::AbstractCamera,
emitter::AbstractEmitter;
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF intensity over camera pixels with optional support region optimization.
This version takes a camera object instead of explicit pixel edges.
Arguments
psf::AbstractPSF
: Point spread function to integratecamera::AbstractCamera
: Camera geometry defining pixel edges in micronsemitter::AbstractEmitter
: Emitter with position in microns relative to camerasupport
: Region to calculate (default: Inf = full image)- If Real: radius in microns around emitter
- If Tuple: explicit (xmin, xmax, ymin, ymax) in microns
sampling::Integer=2
: Number of samples per pixel in each dimensionthreaded::Bool=true
: Whether to use multi-threading for integration
Returns
- Array of pixel values with dimensions y, x
- Values are normalized to sum to 1
- Array indices start at [1,1] for top-left pixel
See also: integrate_pixels_amplitude
, AbstractPSF
integrate_pixels(
psf::AbstractPSF,
pixel_edges_x::AbstractVector,
pixel_edges_y::AbstractVector,
emitters::Vector{<:AbstractEmitter};
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF intensity over camera pixels for multiple emitters with optional support region optimization.
Arguments
psf
: Point spread function to integratepixel_edges_x
,pixel_edges_y
: Arrays defining pixel edge coordinatesemitters
: Vector of emitters with position informationsupport
: Region to calculate for each emitter (default: Inf = full image)- If Real: radius in microns around each emitter
- If Tuple: explicit (xmin, xmax, ymin, ymax) in microns (fixed region for all emitters)
sampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true) Set to false when using with automatic differentiation frameworks
Returns
- Array of integrated PSF intensities with dimensions matching the camera
- Values represent actual photon counts from all emitters
Notes
- Results are the sum of individual emitter contributions (incoherent addition)
- For coherent addition, use
integrate_pixels_amplitude
and sum complex amplitudes
integrate_pixels(
psf::AbstractPSF,
camera::AbstractCamera,
emitters::Vector{<:AbstractEmitter};
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF intensity over camera pixels for multiple emitters with optional support region optimization. This version takes a camera object instead of explicit pixel edges.
Arguments
psf
: Point spread function to integratecamera
: Camera geometry defining pixel edgesemitters
: Vector of emitters with position informationsupport
: Region to calculate for each emitter (default: Inf = full image)sampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
- Array of integrated PSF intensities with dimensions matching the camera
- Values represent actual photon counts from all emitters
integrate_pixels(
psf::SplinePSF,
camera::AbstractCamera,
emitter::AbstractEmitter;
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF over camera pixels using interpolation.
Arguments
psf
: SplinePSF instancecamera
: Camera geometryemitter
: Emitter with position informationsupport
: Region to calculate (default: Inf = full image)- If Real: radius in microns around emitter
- If Tuple: explicit (xmin, xmax, ymin, ymax) in microns
sampling
: Subpixel sampling density for integration accuracythreaded
: Whether to use multi-threading for integration (default: true)
Returns
- Array of integrated PSF intensities with dimensions [ny, nx]
- Values represent actual photon counts based on emitter's photon value
Notes
- For 3D SplinePSFs (when z_range is defined), requires an emitter with a z-coordinate
MicroscopePSFs.integrate_pixels_amplitude
— Functionintegrate_pixels_amplitude(
psf::AbstractPSF,
pixel_edges_x::AbstractVector,
pixel_edges_y::AbstractVector,
emitter::AbstractEmitter;
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF complex amplitude over camera pixels with optional support region optimization.
For each pixel in the specified region, numerically integrates the PSF amplitude using the specified sampling density. Unlike intensity integration, returns unnormalized complex amplitudes which can be used for coherent calculations. Automatically handles z-coordinate if both PSF and emitter support it.
Arguments
psf::AbstractPSF
: Point spread function to integratepixel_edges_x::AbstractVector
: X-coordinate edges of pixels in micronspixel_edges_y::AbstractVector
: Y-coordinate edges of pixels in micronsemitter::AbstractEmitter
: Emitter with position in microns relative to camerasupport
: Region to calculate (default: Inf = full image)- If Real: radius in microns around emitter
- If Tuple: explicit (xmin, xmax, ymin, ymax) in microns
sampling::Integer=2
: Number of samples per pixel in each dimensionthreaded::Bool=true
: Whether to use multi-threading for integration
Returns
Matrix{Complex{T}}
: Integrated complex amplitudes where T matches emitter.photons type- Array is indexed as [y,x] with [1,1] at top-left pixel
- Values are not normalized to preserve complex amplitude relationships
Notes
- For coherent calculations, use this function instead of
integrate_pixels
- Return type is complex to support PSFs with phase information
- To get intensity from amplitude:
abs2.(integrate_pixels_amplitude(...))
See also: integrate_pixels
, AbstractPSF
integrate_pixels_amplitude(
psf::AbstractPSF,
camera::AbstractCamera,
emitter::AbstractEmitter;
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF complex amplitude over camera pixels with optional support region optimization.
This version takes a camera object instead of explicit pixel edges.
Arguments
psf::AbstractPSF
: Point spread function to integratecamera::AbstractCamera
: Camera geometry defining pixel edges in micronsemitter::AbstractEmitter
: Emitter with position in microns relative to camerasupport
: Region to calculate (default: Inf = full image)- If Real: radius in microns around emitter
- If Tuple: explicit (xmin, xmax, ymin, ymax) in microns
sampling::Integer=2
: Number of samples per pixel in each dimensionthreaded::Bool=true
: Whether to use multi-threading for integration
Returns
- Matrix of complex amplitudes
- Array indices start at [1,1] for top-left pixel
See also: integrate_pixels
, AbstractPSF
integrate_pixels_amplitude(
psf::AbstractPSF,
pixel_edges_x::AbstractVector,
pixel_edges_y::AbstractVector,
emitters::Vector{<:AbstractEmitter};
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF complex amplitude over camera pixels for multiple emitters with optional support region optimization.
Arguments
psf
: Point spread function to integratepixel_edges_x
,pixel_edges_y
: Arrays defining pixel edge coordinatesemitters
: Vector of emitters with position informationsupport
: Region to calculate for each emitter (default: Inf = full image)sampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
- Array of integrated complex amplitudes
- Values represent coherently summed field contributions
Notes
- Results are the coherent sum of individual emitter field contributions
- For incoherent addition, use
abs2.(integrate_pixels_amplitude(...))
or useintegrate_pixels
integrate_pixels_amplitude(
psf::AbstractPSF,
camera::AbstractCamera,
emitters::Vector{<:AbstractEmitter};
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF complex amplitude over camera pixels for multiple emitters with optional support region optimization. This version takes a camera object instead of explicit pixel edges.
Arguments
psf
: Point spread function to integratecamera
: Camera geometry defining pixel edgesemitters
: Vector of emitters with position informationsupport
: Region to calculate for each emitter (default: Inf = full image)sampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
- Array of integrated complex amplitudes
- Values represent coherently summed field contributions
integrate_pixels_amplitude(
psf::VectorPSF,
camera::AbstractCamera,
emitter::AbstractEmitter;
sampling::Integer=2,
threaded::Bool=true
)
Specialized version of amplitude integration for VectorPSF that preserves polarization components.
Arguments
psf
: VectorPSF instancecamera
: Camera geometry defining pixel edgesemitter
: Emitter with position informationsampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
- 3D array of integrated complex amplitudes with dimensions [y, x, pol] where pol index 1 = Ex and pol index 2 = Ey
integrate_pixels_amplitude(
psf::VectorPSF,
camera::AbstractCamera,
emitters::Vector{<:AbstractEmitter};
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF complex amplitude over camera pixels for multiple emitters. Special version for VectorPSF that preserves polarization components.
Arguments
psf
: VectorPSF instancecamera
: Camera geometry defining pixel edgesemitters
: Vector of emitters with position informationsampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
- 3D array of integrated complex amplitudes with dimensions [y, x, pol] where pol index 1 = Ex and pol index 2 = Ey
Notes
- Results are the coherent sum of individual emitter field contributions
integrate_pixels_amplitude(
psf::SplinePSF,
camera::AbstractCamera,
emitter::AbstractEmitter;
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF amplitude (complex) over camera pixels.
Arguments
psf
: SplinePSF instancecamera
: Camera geometryemitter
: Emitter with position informationsupport
: Region to calculate (default: Inf = full image)- If Real: radius in microns around emitter
- If Tuple: explicit (xmin, xmax, ymin, ymax) in microns
sampling
: Subpixel sampling density for integration accuracythreaded
: Whether to use multi-threading for integration (default: true)
Returns
- Array of integrated PSF complex amplitudes with dimensions [ny, nx]
Notes
- For 3D SplinePSFs (when z_range is defined), requires an emitter with a z-coordinate
Multi-Emitter Integration
The following functions support integrating PSFs for multiple emitters:
MicroscopePSFs.integrate_pixels
— Methodintegrate_pixels(
psf::AbstractPSF,
camera::AbstractCamera,
emitters::Vector{<:AbstractEmitter};
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF intensity over camera pixels for multiple emitters with optional support region optimization. This version takes a camera object instead of explicit pixel edges.
Arguments
psf
: Point spread function to integratecamera
: Camera geometry defining pixel edgesemitters
: Vector of emitters with position informationsupport
: Region to calculate for each emitter (default: Inf = full image)sampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
- Array of integrated PSF intensities with dimensions matching the camera
- Values represent actual photon counts from all emitters
MicroscopePSFs.integrate_pixels_amplitude
— Methodintegrate_pixels_amplitude(
psf::AbstractPSF,
camera::AbstractCamera,
emitters::Vector{<:AbstractEmitter};
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}} = Inf,
sampling::Integer=2,
threaded::Bool=true
)
Integrate PSF complex amplitude over camera pixels for multiple emitters with optional support region optimization. This version takes a camera object instead of explicit pixel edges.
Arguments
psf
: Point spread function to integratecamera
: Camera geometry defining pixel edgesemitters
: Vector of emitters with position informationsupport
: Region to calculate for each emitter (default: Inf = full image)sampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
- Array of integrated complex amplitudes
- Values represent coherently summed field contributions
Pupil Functions
MicroscopePSFs.PupilFunction
— TypePupilFunction{T}
Represents a pupil function with physical parameters.
Fields
nₐ
: Numerical apertureλ
: Wavelength in μmn
: Refractive indexfield
: Complex-valued pupil function array
MicroscopePSFs.VectorPupilFunction
— TypeVectorPupilFunction{T}
Vector pupil function storing Ex,Ey field components as PupilFunctions.
Fields
nₐ::T
: Numerical apertureλ::T
: Wavelength in μmn_medium::T
: Refractive index of the sample mediumn_coverslip::T
: Refractive index of the coverslipn_immersion::T
: Refractive index of the immersion mediumEx::PupilFunction{T}
: x-component of electric fieldEy::PupilFunction{T}
: y-component of electric field
Zernike Module
MicroscopePSFs.Zernike.ZernikeCoefficients
— TypeZernikeCoefficients{T<:Real}
Mutable structure to hold Zernike coefficients for both magnitude and phase of a pupil function. Uses Noll indexing convention starting at index 1.
Fields
mag::Vector{T}
: Coefficients for magnitude (1-indexed per Noll convention)phase::Vector{T}
: Coefficients for phase (1-indexed per Noll convention)
Notes
- First coefficient (index 1) typically represents piston
- Magnitude coefficients are typically normalized with mag[1] = 1
- Phase coefficients represent phase in radians
- RMS normalization (Noll convention) is used so coefficients directly represent RMS wavefront error
Zernike Polynomial Functions
MicroscopePSFs.Zernike.zernikepolynomial
— Functionzernikepolynomial(n::Integer, l::Integer, ρ::Real, ϕ::Real) -> Real
Compute the complete Zernike polynomial Z_n^l(ρ,ϕ) with Noll normalization (RMS=1).
Arguments
n
: Radial orderl
: Azimuthal frequency (signed)ρ
: Radial coordinate (0 ≤ ρ ≤ 1)ϕ
: Azimuthal angle in radians
Notes
- Uses Noll normalization where RMS=1 over unit circle
- Combines radial polynomial with appropriate trigonometric function
zernikepolynomial(j::Integer, ρ::Real, ϕ::Real) -> Real
Compute Zernike polynomial using Noll index.
Arguments
j
: Polynomial index using Noll conventionρ
: Radial coordinateϕ
: Azimuthal angle
MicroscopePSFs.Zernike.radialpolynomial
— Functionradialpolynomial(n::Integer, m::Integer, ρ::Real) -> Real
Compute the unnormalized radial component R_n^m(ρ) of the Zernike polynomial.
Arguments
n
: Radial orderm
: Azimuthal order (absolute value of azimuthal frequency)ρ
: Radial coordinate (0 ≤ ρ ≤ 1)
Notes
- Returns 0 for ρ > 1
- No normalization applied - returns the standard mathematical form
MicroscopePSFs.Zernike.max_radial_order
— Functionmax_radial_order(num_coeffs::Integer)
Calculate maximum radial order N given number of coefficients L. Solves (N+1)(N+2)/2 = L
MicroscopePSFs.Zernike.evaluate_pupil
— Functionevaluate_pupil(coeffs::ZernikeCoefficients, grid_size::Integer) -> Matrix{Complex{Float64}}
Generate complex pupil function from Zernike coefficients.
Arguments
coeffs
: ZernikeCoefficients containing magnitude and phase coefficientsgrid_size
: Size of the output grid (gridsize × gridsize)
Returns
- Complex-valued matrix representing the pupil function
Notes
- Output grid is normalized to unit circle
- Points outside unit circle are set to zero
- Phase is applied as exp(iϕ)
- Uses Noll indexing throughout
Zernike Analysis Functions
MicroscopePSFs.Zernike.rms
— Functionrms(coeffs::ZernikeCoefficients) -> Tuple{Float64,Float64}
Calculate RMS values for magnitude and phase coefficients (excluding piston).
MicroscopePSFs.Zernike.significant_terms
— Functionsignificant_terms(coeffs::ZernikeCoefficients,
threshold::Real=0.01) -> Vector{Tuple{Int,Float64,Float64}}
Return list of significant terms: (index, magnitude, phase) above threshold.
Index Conversion
MicroscopePSFs.Zernike.nl2osa
— Functionnl2osa(n::Integer, l::Integer) -> Int
Convert from (n,l) indices to OSA/ANSI single index.
Arguments
n
: Radial order (≥ 0)l
: Azimuthal frequency, must satisfy:- |l| ≤ n
- n - |l| must be even
Returns
- OSA/ANSI single index j = (n(n+2) + l)/2
Examples
julia> nl2osa(4, 2)
15
MicroscopePSFs.Zernike.osa2nl
— Functionosa2nl(j::Integer) -> Tuple{Int,Int}
Convert from OSA/ANSI single index to (n,l) indices.
Arguments
j
: OSA/ANSI index (≥ 0)
Returns
- Tuple of (n,l) indices
MicroscopePSFs.Zernike.nl2noll
— Functionnl2noll(n::Integer, l::Integer) -> Int
Convert from (n,l) indices to Noll single index.
Arguments
n
: Radial order (≥ 0)l
: Azimuthal frequency, must satisfy:- |l| ≤ n
- n - |l| must be even
Returns
- Noll single index
MicroscopePSFs.Zernike.noll2nl
— Functionnoll2nl(j::Integer) -> Tuple{Int,Int}
Convert from Noll single index to (n,l) indices.
Arguments
j
: Noll index (> 0)
Returns
- Tuple of (n,l) indices
MicroscopePSFs.Zernike.osa2noll
— Functionosa2noll(j::Integer) -> Int
Convert from OSA/ANSI to Noll index.
MicroscopePSFs.Zernike.noll2osa
— Functionnoll2osa(j::Integer) -> Int
Convert from Noll to OSA/ANSI index.
Emitters
MicroscopePSFs.DipoleVector
— TypeDipoleVector{T} <: Real
A 3D vector representing the dipole orientation.
Fields
px::T
: x component of the dipole vectorpy::T
: y component of the dipole vectorpz::T
: z component of the dipole vector
MicroscopePSFs.DipoleEmitter3D
— TypeDipoleEmitter3D{T} <: AbstractEmitter
3D dipole emitter with position, orientation and optical properties.
Fields
x::T
: x-coordinate in micronsy::T
: y-coordinate in micronsz::T
: z-coordinate in micronsphotons::T
: number of photonsdipole::DipoleVector{T}
: dipole orientation vector
I/O Functions
MicroscopePSFs.save_psf
— Functionsave_psf(filename::String, object; metadata::Dict=Dict())
Save a PSF or related object (e.g., ZernikeCoefficients, PupilFunction) to an HDF5 file.
Arguments
filename
: Path where the PSF will be savedobject
: Object to save (PSF, ZernikeCoefficients, PupilFunction, etc.)metadata
: Optional dictionary of additional metadata to include
Returns
filename
for chaining
MicroscopePSFs.load_psf
— Functionload_psf(filename::String)
Load a PSF or related object (e.g., ZernikeCoefficients, PupilFunction) from an HDF5 file.
Arguments
filename
: Path to the HDF5 file containing a saved PSF object
Returns
- The loaded object with its original type (PSF, ZernikeCoefficients, PupilFunction, etc.)
Examples
# Load a previously saved PSF
psf = load_psf("my_airy_psf.h5")
# Use the loaded PSF normally
intensity = psf(0.1, 0.2)
# Load Zernike coefficients and use them
zc = load_psf("zernike_coeffs.h5")
psf = ScalarPSF(1.4, 0.532, 1.518; zernike_coeffs=zc)
See also: save_psf
Complete API (All Documented Functions)
This section lists additional internal functions and types that are documented but not part of the public API.
MicroscopePSFs.Abstract2DPSF
— TypeAbstract2DPSF{T<:AbstractFloat}
Abstract type for all 2D point spread functions. Parameterized by numeric type T.
MicroscopePSFs.Abstract3DPSF
— TypeAbstract3DPSF{T<:AbstractFloat}
Abstract type for all 3D point spread functions.
MicroscopePSFs._calculate_field_amplitude
— Method_calculate_field_amplitude(
pupil::VectorPupilFunction,
nₐ::Real, λ::Real, n_medium::Real, n_coverslip::Real, n_immersion::Real, z_stage::Real,
x::Real, y::Real, z::Real)
Internal helper function to compute complex vector amplitude at given position.
Arguments
pupil
: Vector pupil function containing Ex and Ey fieldsnₐ
,λ
,n_medium
,n_coverslip
,n_immersion
: Optical parametersz_stage
: Distance the sample stage was moved away from the nominal focal planex, y, z
: Position in microns
Returns
- Vector [Ex, Ey] of complex field amplitudes
Notes
- This is a helper function used by both amplitude() and the PSF evaluation function
- Implements the core field calculation for a single pupil
MicroscopePSFs._check_normalization
— Method_check_normalization(values; tol=1e-6)
Internal function to verify array sums to 1 within tolerance.
MicroscopePSFs._integrate_pixels_core!
— Method_integrate_pixels_core!(
result::AbstractArray,
pixel_edges_x::AbstractVector,
pixel_edges_y::AbstractVector,
func::Function,
emitter_x::Real,
emitter_y::Real;
sampling::Integer=2,
threaded::Bool=true
)
Internal function implementing the core pixel integration logic. Supports both scalar and vector/tensor return types from the function.
Arguments
result
: Pre-allocated array where integration results will be storedpixel_edges_x
,pixel_edges_y
: Arrays defining pixel edge coordinatesfunc
: Function to evaluate at PSF coordinates - can return scalar, vector, or tensoremitter_x
,emitter_y
: Emitter coordinates in same units as pixel edgessampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true) Set to false when using with automatic differentiation frameworks
Returns
result
array filled with integration values
MicroscopePSFs._integrate_pixels_generic!
— Method_integrate_pixels_generic!(
result::AbstractArray,
psf::AbstractPSF,
pixel_edges_x::AbstractVector,
pixel_edges_y::AbstractVector,
emitter::AbstractEmitter,
f::Function;
sampling::Integer=2,
threaded::Bool=true
)
Generic wrapper for PSF integration that handles coordinate systems and validation. Automatically handles z-coordinates when both PSF and emitter support them. Supports functions that return either scalars or vectors/tensors.
Arguments
result
: Pre-allocated array where results will be storedpsf
: Point spread function to integratepixel_edges_x
,pixel_edges_y
: Arrays defining pixel edge coordinatesemitter
: Emitter with position informationf
: Function to evaluate (e.g., (p,x,y) -> p(x,y) for intensity)sampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
result
array filled with integration values
MicroscopePSFs._integrate_pixels_generic
— Method_integrate_pixels_generic(
psf::AbstractPSF,
pixel_edges_x::AbstractVector,
pixel_edges_y::AbstractVector,
emitter::AbstractEmitter,
f::Function,
::Type{T},
output_dims=();
sampling::Integer=2,
threaded::Bool=true
) where T
Internal generic integration routine used by both intensity and amplitude integration. Creates a new result array and delegates to the in-place version. Supports both scalar and vector/tensor return types.
Arguments
output_dims
: Additional dimensions for the result array beyond the basic [y,x] dimensions. For vector returns (e.g., [Ex, Ey]), set output_dims=(2,) for a [y,x,2] result array.threaded
: Whether to use multi-threading for integration (default: true)
Returns
- Array with dimensions [y,x] or [y,x,output_dims...] depending on the function's return type
MicroscopePSFs._io_check_required_fields
— Method_io_check_required_fields(group, fields)
Check if all required fields exist in an HDF5 group. Throws an error if any required field is missing.
MicroscopePSFs._io_load_complex_array
— Method_io_load_complex_array(group, name)
Load a complex array from an HDF5 group by combining real and imaginary parts.
MicroscopePSFs._io_load_range
— Method_io_load_range(group, name)
Load a range from an HDF5 group by reconstructing from start, step, and length.
MicroscopePSFs._io_load_zernike_coeffs
— Method_io_load_zernike_coeffs(file; group_name="zernike_coefficients")
Load ZernikeCoefficients from an HDF5 file from a specified group. Returns nothing if the group doesn't exist.
MicroscopePSFs._io_save_complex_array
— Method_io_save_complex_array(group, name, array)
Save a complex array to an HDF5 group by splitting into real and imaginary parts.
MicroscopePSFs._io_save_pupil_params
— Method_io_save_pupil_params(params, pupil)
Save common pupil parameters to an HDF5 group.
MicroscopePSFs._io_save_range
— Method_io_save_range(group, name, range)
Save a range to an HDF5 group by storing start, step, and length.
MicroscopePSFs._io_save_zernike_coeffs
— Method_io_save_zernike_coeffs(file, coeffs; group_name="zernike_coefficients")
Save ZernikeCoefficients to an HDF5 file in a specified group. Does nothing if coeffs is nothing.
MicroscopePSFs._load_psf_impl
— Method_load_psf_impl(file::HDF5.File, ::Type{AiryPSF})
Load an AiryPSF from an HDF5 file.
MicroscopePSFs._load_psf_impl
— Method_load_psf_impl(file::HDF5.File, ::Type{GaussianPSF})
Load a GaussianPSF from an HDF5 file.
MicroscopePSFs._load_psf_impl
— Method_load_psf_impl(file::HDF5.File, ::Type{PupilFunction})
Load a PupilFunction from an HDF5 file.
MicroscopePSFs._load_psf_impl
— Method_load_psf_impl(file::HDF5.File, ::Type{ScalarPSF})
Load a ScalarPSF from an HDF5 file.
MicroscopePSFs._load_psf_impl
— Method_load_psf_impl(file::HDF5.File, ::Type{SplinePSF})
Load a SplinePSF from an HDF5 file, reconstructing it using the standard constructor.
MicroscopePSFs._load_psf_impl
— Method_load_psf_impl(file::HDF5.File, ::Type{VectorPSF})
Load a VectorPSF from an HDF5 file, reconstructing all components.
Arguments
file
: Open HDF5 file handle::Type{VectorPSF}
: Type to load
Returns
- VectorPSF reconstructed from stored data
MicroscopePSFs._load_psf_impl
— Method_load_psf_impl(file::HDF5.File, ::Type{VectorPupilFunction})
Load a VectorPupilFunction from an HDF5 file.
MicroscopePSFs._load_psf_impl
— Method_load_psf_impl(file::HDF5.File, ::Type{ZernikeCoefficients})
Load ZernikeCoefficients from an HDF5 file.
MicroscopePSFs._sample_psf_2d
— Method_sample_psf_2d(psf::AbstractPSF, x_range::AbstractRange, y_range::AbstractRange)
Internal helper function to sample a 2D PSF on a regular grid. Returns an array with dimensions [y, x] containing raw PSF values.
MicroscopePSFs._sample_psf_3d
— Method_sample_psf_3d(psf::AbstractPSF, x_range::AbstractRange, y_range::AbstractRange, z_range::AbstractRange)
Internal helper function to sample a 3D PSF on a regular grid. Returns an array with dimensions [y, x, z] containing raw PSF values.
MicroscopePSFs._save_psf_impl
— Method_save_psf_impl(file::HDF5.File, psf::AiryPSF)
Save an AiryPSF to an HDF5 file.
MicroscopePSFs._save_psf_impl
— Method_save_psf_impl(file::HDF5.File, psf::GaussianPSF)
Save a GaussianPSF to an HDF5 file.
MicroscopePSFs._save_psf_impl
— Method_save_psf_impl(file::HDF5.File, pupil::PupilFunction)
Save a PupilFunction to an HDF5 file, including physical parameters and complex field.
MicroscopePSFs._save_psf_impl
— Method_save_psf_impl(file::HDF5.File, psf::ScalarPSF)
Save a ScalarPSF to an HDF5 file, including pupil function and any Zernike coefficients.
MicroscopePSFs._save_psf_impl
— Method_save_psf_impl(file::HDF5.File, psf::SplinePSF)
Save a SplinePSF to an HDF5 file, including the original grid data used to create the spline interpolation.
MicroscopePSFs._save_psf_impl
— Method_save_psf_impl(file::HDF5.File, psf::VectorPSF)
Save a VectorPSF to an HDF5 file, including all optical parameters, dipole orientation, and vector pupil fields.
Arguments
file
: Open HDF5 file handlepsf
: VectorPSF to save
MicroscopePSFs._save_psf_impl
— Method_save_psf_impl(file::HDF5.File, pupil::VectorPupilFunction)
Save a VectorPupilFunction to an HDF5 file, including all physical parameters and both Ex and Ey field components.
MicroscopePSFs._save_psf_impl
— Method_save_psf_impl(file::HDF5.File, zc::ZernikeCoefficients)
Save ZernikeCoefficients to an HDF5 file, storing magnitude and phase arrays.
MicroscopePSFs.apply_aperture!
— Functionapply_aperture!(p::PupilFunction, radius::Real=1.0)
Apply circular aperture to pupil function. Radius is relative to NA.
MicroscopePSFs.apply_defocus!
— Methodapply_defocus!(p::PupilFunction, z::Real)
Apply defocus phase to pupil function for propagation distance z.
MicroscopePSFs.calculate_apodization
— Methodcalculate_apodization(kr2::Real, λ::Real,
n_medium::Real, n_immersion::Real)
Calculate apodization factor for energy conservation.
Arguments
kr2
: Squared lateral spatial frequencyλ
: Wavelength in micronsn_medium
: Refractive index of the sample mediumn_immersion
: Refractive index of the immersion medium
Returns
- Apodization factor for energy conservation
MicroscopePSFs.calculate_axial_phase
— Functioncalculate_axial_phase(z::Real, z_stage::Real, kz_medium::Complex,
kz_coverslip::Complex, kz_immersion::Complex,
coverslip_thickness::Real=0.17)
Calculate total axial phase from defocus.
Arguments
z
: Emitter position relative to the coverslip (depth above the coverslip in μm)z_stage
: Distance the sample stage was moved away from the nominal focal plane at the coverslip (μm)kz_medium
: z-component of wave vector in sample mediumkz_coverslip
: z-component of wave vector in coverslipkz_immersion
: z-component of wave vector in immersion mediumcoverslip_thickness
: Thickness of coverslip in mm (default: 0.17mm)
Returns
- Defocus phase
MicroscopePSFs.calculate_dipole_field_components
— Methodcalculate_dipole_field_components(ϕ::Complex, sinθ::Complex, cosθ::Complex,
dipole::DipoleVector, Tp::Complex, Ts::Complex)
Calculate vectorial field components for a dipole orientation.
Arguments
ϕ
: Azimuthal angle in pupil planesinθ
: Sine of polar anglecosθ
: Cosine of polar angledipole
: Dipole orientation vectorTp
: p-polarization transmission coefficientTs
: s-polarization transmission coefficient
Returns
- Tuple (Ex, Ey) of complex field components
MicroscopePSFs.calculate_interface_fresnel
— Methodcalculate_interface_fresnel(kr2::Real, λ::Real,
n1::Real, n2::Real)
Calculate Fresnel transmission coefficients for a single interface.
Arguments
kr2
: Squared lateral spatial frequencyλ
: Wavelength in micronsn1
: Refractive index of first mediumn2
: Refractive index of second medium
Returns
- Tuple of (Tp, Ts): p and s polarization transmission coefficients
MicroscopePSFs.calculate_wave_vectors
— Methodcalculate_wave_vectors(kr2::Real, λ::Real,
n_medium::Real, n_coverslip::Real, n_immersion::Real)
Calculate z-components of wave vectors in all three media.
Arguments
kr2
: Squared lateral spatial frequencyλ
: Wavelength in micronsn_medium
: Refractive index of the sample mediumn_coverslip
: Refractive index of the coverslipn_immersion
: Refractive index of the immersion medium
Returns
- Tuple of (kzmedium, kzcoverslip, kz_immersion): z-components of wave vectors in each medium
MicroscopePSFs.fill_vector_pupils!
— Functionfill_vector_pupils!(vpupil::VectorPupilFunction, dipole::DipoleVector,
base_pupil::Union{Nothing, PupilFunction}=nothing;
normalize::Bool=true)
Fill vector pupil function with field components including dipole orientation, base aberrations, and proper apodization.
This pre-calculates all position-independent factors of the pupil function.
Arguments
vpupil
: Vector pupil function to filldipole
: Dipole orientation vectorbase_pupil
: Optional base aberration pupil function
Keyword Arguments
normalize
: Whether to normalize the pupil function (default: true)
Returns
- Filled vector pupil function (normalized if
normalize=true
)
MicroscopePSFs.get_pixel_indices
— Methodget_pixel_indices(
pixel_edges_x::AbstractVector,
pixel_edges_y::AbstractVector,
emitter::AbstractEmitter,
support::Union{Real,Tuple{<:Real,<:Real,<:Real,<:Real}}
) -> Tuple{UnitRange{Int},UnitRange{Int}}
Get the pixel indices that overlap with the support region.
Arguments
pixel_edges_x
,pixel_edges_y
: Arrays defining pixel edge coordinatesemitter
: Emitter with position informationsupport
: Region of interest, either a radius or explicit bounds (xmin, xmax, ymin, ymax)
Returns
- Tuple of (irange, jrange) with pixel indices that cover the support region If the support region doesn't overlap with the camera, returns empty ranges
MicroscopePSFs.integrate_pixels!
— Methodintegrate_pixels!(
result::AbstractMatrix{T},
psf::AbstractPSF,
camera::AbstractCamera,
emitter::AbstractEmitter;
sampling::Integer=2,
threaded::Bool=true
) where T <: Real
Integrate PSF intensity over camera pixels using camera object, storing the result in a pre-allocated matrix.
Arguments
result
: Pre-allocated array where results will be storedpsf
: Point spread function to integratecamera
: Camera geometry defining pixel edgesemitter
: Emitter with position informationsampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
- The
result
array, now filled with integrated intensities
MicroscopePSFs.integrate_pixels!
— Methodintegrate_pixels!(
result::AbstractMatrix{T},
psf::AbstractPSF,
pixel_edges_x::AbstractVector,
pixel_edges_y::AbstractVector,
emitter::AbstractEmitter;
sampling::Integer=2,
threaded::Bool=true
) where T <: Real
Integrate PSF intensity over camera pixels, storing the result in a pre-allocated matrix. Automatically uses z-coordinate if both PSF and emitter support it.
Arguments
result
: Pre-allocated array where results will be storedpsf
: Point spread function to integratepixel_edges_x
,pixel_edges_y
: Arrays defining pixel edge coordinatesemitter
: Emitter with position informationsampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true) Set to false when using with automatic differentiation frameworks
Returns
- The
result
array, now filled with integrated intensities scaled by emitter.photons
Notes
- Array is indexed as [y,x] with [1,1] at top-left pixel
MicroscopePSFs.integrate_pixels_amplitude!
— Methodintegrate_pixels_amplitude!(
result::AbstractMatrix{Complex{T}},
psf::AbstractPSF,
camera::AbstractCamera,
emitter::AbstractEmitter;
sampling::Integer=2,
threaded::Bool=true
) where T <: Real
Integrate PSF complex amplitude over camera pixels using camera object, storing the result in a pre-allocated matrix.
Arguments
result
: Pre-allocated complex array where results will be storedpsf
: Point spread function to integratecamera
: Camera geometry defining pixel edgesemitter
: Emitter with position informationsampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
- The
result
array, now filled with integrated complex amplitudes
MicroscopePSFs.integrate_pixels_amplitude!
— Methodintegrate_pixels_amplitude!(
result::AbstractMatrix{Complex{T}},
psf::AbstractPSF,
pixel_edges_x::AbstractVector,
pixel_edges_y::AbstractVector,
emitter::AbstractEmitter;
sampling::Integer=2,
threaded::Bool=true
) where T <: Real
Integrate PSF complex amplitude over camera pixels, storing the result in a pre-allocated matrix. Automatically uses z-coordinate if both PSF and emitter support it.
Arguments
result
: Pre-allocated complex array where results will be storedpsf
: Point spread function to integratepixel_edges_x
,pixel_edges_y
: Arrays defining pixel edge coordinatesemitter
: Emitter with position informationsampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
- The
result
array, now filled with integrated complex amplitudes
Notes
- Array is indexed as [y,x] with [1,1] at top-left pixel
MicroscopePSFs.integrate_pixels_amplitude!
— Methodintegrate_pixels_amplitude!(
result::AbstractArray{Complex{T},3},
psf::VectorPSF,
camera::AbstractCamera,
emitter::AbstractEmitter;
sampling::Integer=2,
threaded::Bool=true
) where T <: Real
Special version of amplitude integration for VectorPSF that preserves polarization components. Uses the flexible core integration function that handles vector returns.
Arguments
result
: Pre-allocated 3D complex array where results will be stored, dimensions [y, x, pol]psf
: VectorPSF instancecamera
: Camera geometry defining pixel edgesemitter
: Emitter with position informationsampling
: Subpixel sampling density (default: 2)threaded
: Whether to use multi-threading for integration (default: true)
Returns
- The
result
array, filled with integrated complex amplitudes for each polarization
MicroscopePSFs.kmax
— Methodkmax(p::PupilFunction)
Get maximum spatial frequency in μm⁻¹
MicroscopePSFs.kmax
— MethodGet maximum spatial frequency in μm⁻¹
MicroscopePSFs.kpixelsize
— Methodkpixelsize(p::PupilFunction)
Get pupil plane sampling in μm⁻¹
MicroscopePSFs.kpixelsize
— MethodGet pupil plane sampling in μm⁻¹
MicroscopePSFs.k₀
— Methodk₀(p::PupilFunction)
Get central wavevector magnitude in μm⁻¹
MicroscopePSFs.k₀
— MethodGet central wavevector magnitude in μm⁻¹
MicroscopePSFs.normalize!
— Methodnormalize!(p::PupilFunction)
Normalize pupil function to unit energy using Parseval's theorem.
MicroscopePSFs.normalize!
— Methodnormalize!(p::VectorPupilFunction)
Normalize the electric field components of the vector pupil function. Ensures total energy across both components equals 1.
Arguments
p
: Vector pupil function to normalize
Returns
- Normalized vector pupil function
MicroscopePSFs.supports_3d
— Methodsupports_3d(psf_type::Type{<:AbstractPSF}) -> Bool
Check if a PSF type supports 3D evaluation via psf(x, y, z).
MicroscopePSFs.update_pupils!
— Methodupdate_pupils!(psf::VectorPSF) -> VectorPSF
Update the vector pupil functions based on stored Zernike coefficients and/or base pupil. This is useful after modifying aberrations to regenerate the pupil fields.
Arguments
psf
: VectorPSF to update
Returns
- Updated VectorPSF
Notes
- Requires either stored Zernike coefficients or a base pupil
- For rotating dipoles (multiple pupils), updates each pupil with the appropriate dipole orientation
- Returns the updated PSF for method chaining
Zernike Module Internal API
MicroscopePSFs.Zernike.get_nl
— Methodget_nl(j::Integer) -> Tuple{Int,Int}
Get (n,l) indices from Noll index.