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

2D PSF Models

MicroscopePSFs.GaussianPSFType
GaussianPSF{T<:AbstractFloat} <: Abstract2DPSF{T}

Isotropic 2D Gaussian PSF.

Fields

  • σ: Standard deviation in physical units (typically microns)
source
MicroscopePSFs.AiryPSFType
AiryPSF{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ₐ/λ
source

3D PSF Models

MicroscopePSFs.ScalarPSFType
ScalarPSF{T} <: Abstract3DPSF{T}

Scalar 3D PSF using explicit pupil function representation.

Fields

  • nₐ::T: Numerical aperture
  • λ::T: Wavelength in microns
  • n::T: Refractive index
  • pupil::PupilFunction{T}: Complex pupil function
  • zernike_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.

source
MicroscopePSFs.VectorPSFType
VectorPSF{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 microns
  • n_medium::T: Sample medium refractive index
  • n_coverslip::T: Cover slip refractive index
  • n_immersion::T: Immersion medium refractive index
  • dipole::DipoleVector{T}: Dipole orientation
  • z_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 aberrations
  • zernike_coeffs::Union{Nothing, ZernikeCoefficients{T}}: Zernike coefficients used to create this PSF (if applicable)
source

PSF Acceleration

MicroscopePSFs.SplinePSFType
SplinePSF{T<:AbstractFloat, IT<:AbstractInterpolation} <: AbstractPSF

A point spread function (PSF) represented as a B-spline interpolation.

Fields

  • spline: The B-spline interpolation object
  • x_range: Range of x-coordinates used for uniform grid interpolation
  • y_range: Range of y-coordinates used for uniform grid interpolation
  • z_range: Range of z-coordinates for 3D PSFs, or nothing for 2D PSFs
  • original_grid: Original grid data used to create the interpolation
  • interp_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
source

Core Interface Functions

MicroscopePSFs.amplitudeFunction
amplitude(p::PupilFunction, x::Real, y::Real, z::Real)

Calculate complex amplitude from pupil function integration.

Arguments

  • p::PupilFunction: Pupil function
  • x::Real: X position in μm
  • y::Real: Y position in μm
  • z::Real: Z position in μm
source
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
source
amplitude(psf::AbstractPSF, x::Real, y::Real, z::Real)

Evaluate complex field amplitude in 3D. z is axial distance from focus in microns.

source
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 instance
  • x, 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
source
amplitude(psf::VectorPSF, x::Real, y::Real, z::Real)

Compute complex vector amplitude at given position.

Arguments

  • psf: Vector PSF instance
  • x, y: Lateral position in microns relative to PSF center
  • z: 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)
source
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 instance
  • x, 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
source
amplitude(psf::SplinePSF, x::Real, y::Real)

Calculate the complex amplitude of the PSF at position (x, y) with z=0.

Arguments

  • psf: SplinePSF instance
  • x, y: Coordinates in microns relative to PSF center

Returns

  • Complex amplitude = sqrt(intensity) with zero phase
source
MicroscopePSFs.integrate_pixelsFunction
integrate_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 integrate
  • pixel_edges_x::AbstractVector: X-coordinate edges of pixels in microns
  • pixel_edges_y::AbstractVector: Y-coordinate edges of pixels in microns
  • emitter::AbstractEmitter: Emitter with position in microns relative to camera
  • support: 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 dimension
  • threaded::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

source
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 integrate
  • camera::AbstractCamera: Camera geometry defining pixel edges in microns
  • emitter::AbstractEmitter: Emitter with position in microns relative to camera
  • support: 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 dimension
  • threaded::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

source
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 integrate
  • pixel_edges_x, pixel_edges_y: Arrays defining pixel edge coordinates
  • emitters: Vector of emitters with position information
  • support: 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
source
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 integrate
  • camera: Camera geometry defining pixel edges
  • emitters: Vector of emitters with position information
  • support: 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
source
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 instance
  • camera: Camera geometry
  • emitter: Emitter with position information
  • support: 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 accuracy
  • threaded: 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
source
MicroscopePSFs.integrate_pixels_amplitudeFunction
integrate_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 integrate
  • pixel_edges_x::AbstractVector: X-coordinate edges of pixels in microns
  • pixel_edges_y::AbstractVector: Y-coordinate edges of pixels in microns
  • emitter::AbstractEmitter: Emitter with position in microns relative to camera
  • support: 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 dimension
  • threaded::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

source
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 integrate
  • camera::AbstractCamera: Camera geometry defining pixel edges in microns
  • emitter::AbstractEmitter: Emitter with position in microns relative to camera
  • support: 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 dimension
  • threaded::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

source
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 integrate
  • pixel_edges_x, pixel_edges_y: Arrays defining pixel edge coordinates
  • emitters: Vector of emitters with position information
  • support: 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 use integrate_pixels
source
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 integrate
  • camera: Camera geometry defining pixel edges
  • emitters: Vector of emitters with position information
  • support: 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
source
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 instance
  • camera: Camera geometry defining pixel edges
  • emitter: Emitter with position information
  • sampling: 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
source
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 instance
  • camera: Camera geometry defining pixel edges
  • emitters: Vector of emitters with position information
  • sampling: 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
source
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 instance
  • camera: Camera geometry
  • emitter: Emitter with position information
  • support: 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 accuracy
  • threaded: 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
source

Multi-Emitter Integration

The following functions support integrating PSFs for multiple emitters:

MicroscopePSFs.integrate_pixelsMethod
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 integrate
  • camera: Camera geometry defining pixel edges
  • emitters: Vector of emitters with position information
  • support: 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
source
MicroscopePSFs.integrate_pixels_amplitudeMethod
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 integrate
  • camera: Camera geometry defining pixel edges
  • emitters: Vector of emitters with position information
  • support: 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
source

Pupil Functions

MicroscopePSFs.PupilFunctionType
PupilFunction{T}

Represents a pupil function with physical parameters.

Fields

  • nₐ: Numerical aperture
  • λ: Wavelength in μm
  • n: Refractive index
  • field: Complex-valued pupil function array
source
MicroscopePSFs.VectorPupilFunctionType
VectorPupilFunction{T}

Vector pupil function storing Ex,Ey field components as PupilFunctions.

Fields

  • nₐ::T: Numerical aperture
  • λ::T: Wavelength in μm
  • n_medium::T: Refractive index of the sample medium
  • n_coverslip::T: Refractive index of the coverslip
  • n_immersion::T: Refractive index of the immersion medium
  • Ex::PupilFunction{T}: x-component of electric field
  • Ey::PupilFunction{T}: y-component of electric field
source

Zernike Module

MicroscopePSFs.Zernike.ZernikeCoefficientsType
ZernikeCoefficients{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
source

Zernike Polynomial Functions

MicroscopePSFs.Zernike.zernikepolynomialFunction
zernikepolynomial(n::Integer, l::Integer, ρ::Real, ϕ::Real) -> Real

Compute the complete Zernike polynomial Z_n^l(ρ,ϕ) with Noll normalization (RMS=1).

Arguments

  • n: Radial order
  • l: 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
source
zernikepolynomial(j::Integer, ρ::Real, ϕ::Real) -> Real

Compute Zernike polynomial using Noll index.

Arguments

  • j: Polynomial index using Noll convention
  • ρ: Radial coordinate
  • ϕ: Azimuthal angle
source
MicroscopePSFs.Zernike.radialpolynomialFunction
radialpolynomial(n::Integer, m::Integer, ρ::Real) -> Real

Compute the unnormalized radial component R_n^m(ρ) of the Zernike polynomial.

Arguments

  • n: Radial order
  • m: Azimuthal order (absolute value of azimuthal frequency)
  • ρ: Radial coordinate (0 ≤ ρ ≤ 1)

Notes

  • Returns 0 for ρ > 1
  • No normalization applied - returns the standard mathematical form
source
MicroscopePSFs.Zernike.evaluate_pupilFunction
evaluate_pupil(coeffs::ZernikeCoefficients, grid_size::Integer) -> Matrix{Complex{Float64}}

Generate complex pupil function from Zernike coefficients.

Arguments

  • coeffs: ZernikeCoefficients containing magnitude and phase coefficients
  • grid_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
source

Zernike Analysis Functions

MicroscopePSFs.Zernike.rmsFunction
rms(coeffs::ZernikeCoefficients) -> Tuple{Float64,Float64}

Calculate RMS values for magnitude and phase coefficients (excluding piston).

source
MicroscopePSFs.Zernike.significant_termsFunction
significant_terms(coeffs::ZernikeCoefficients, 
                 threshold::Real=0.01) -> Vector{Tuple{Int,Float64,Float64}}

Return list of significant terms: (index, magnitude, phase) above threshold.

source

Index Conversion

MicroscopePSFs.Zernike.nl2osaFunction
nl2osa(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
source
MicroscopePSFs.Zernike.osa2nlFunction
osa2nl(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
source
MicroscopePSFs.Zernike.nl2nollFunction
nl2noll(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
source
MicroscopePSFs.Zernike.noll2nlFunction
noll2nl(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
source

Emitters

MicroscopePSFs.DipoleVectorType
DipoleVector{T} <: Real

A 3D vector representing the dipole orientation.

Fields

  • px::T: x component of the dipole vector
  • py::T: y component of the dipole vector
  • pz::T: z component of the dipole vector
source
MicroscopePSFs.DipoleEmitter3DType
DipoleEmitter3D{T} <: AbstractEmitter

3D dipole emitter with position, orientation and optical properties.

Fields

  • x::T: x-coordinate in microns
  • y::T: y-coordinate in microns
  • z::T: z-coordinate in microns
  • photons::T: number of photons
  • dipole::DipoleVector{T}: dipole orientation vector
source

I/O Functions

MicroscopePSFs.save_psfFunction
save_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 saved
  • object: Object to save (PSF, ZernikeCoefficients, PupilFunction, etc.)
  • metadata: Optional dictionary of additional metadata to include

Returns

  • filename for chaining
source
MicroscopePSFs.load_psfFunction
load_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

source

Complete API (All Documented Functions)

This section lists additional internal functions and types that are documented but not part of the public API.

MicroscopePSFs._calculate_field_amplitudeMethod
_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 fields
  • nₐ, λ, n_medium, n_coverslip, n_immersion: Optical parameters
  • z_stage: Distance the sample stage was moved away from the nominal focal plane
  • x, 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
source
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 stored
  • pixel_edges_x, pixel_edges_y: Arrays defining pixel edge coordinates
  • func: Function to evaluate at PSF coordinates - can return scalar, vector, or tensor
  • emitter_x, emitter_y: Emitter coordinates in same units as pixel edges
  • 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

  • result array filled with integration values
source
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 stored
  • psf: Point spread function to integrate
  • pixel_edges_x, pixel_edges_y: Arrays defining pixel edge coordinates
  • emitter: Emitter with position information
  • f: 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
source
MicroscopePSFs._integrate_pixels_genericMethod
_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
source
MicroscopePSFs._io_load_zernike_coeffsMethod
_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.

source
MicroscopePSFs._io_save_zernike_coeffsMethod
_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.

source
MicroscopePSFs._load_psf_implMethod
_load_psf_impl(file::HDF5.File, ::Type{SplinePSF})

Load a SplinePSF from an HDF5 file, reconstructing it using the standard constructor.

source
MicroscopePSFs._load_psf_implMethod
_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
source
MicroscopePSFs._sample_psf_2dMethod
_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.

source
MicroscopePSFs._sample_psf_3dMethod
_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.

source
MicroscopePSFs._save_psf_implMethod
_save_psf_impl(file::HDF5.File, pupil::PupilFunction)

Save a PupilFunction to an HDF5 file, including physical parameters and complex field.

source
MicroscopePSFs._save_psf_implMethod
_save_psf_impl(file::HDF5.File, psf::ScalarPSF)

Save a ScalarPSF to an HDF5 file, including pupil function and any Zernike coefficients.

source
MicroscopePSFs._save_psf_implMethod
_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.

source
MicroscopePSFs._save_psf_implMethod
_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 handle
  • psf: VectorPSF to save
source
MicroscopePSFs._save_psf_implMethod
_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.

source
MicroscopePSFs._save_psf_implMethod
_save_psf_impl(file::HDF5.File, zc::ZernikeCoefficients)

Save ZernikeCoefficients to an HDF5 file, storing magnitude and phase arrays.

source
MicroscopePSFs.calculate_apodizationMethod
calculate_apodization(kr2::Real, λ::Real, 
                    n_medium::Real, n_immersion::Real)

Calculate apodization factor for energy conservation.

Arguments

  • kr2: Squared lateral spatial frequency
  • λ: Wavelength in microns
  • n_medium: Refractive index of the sample medium
  • n_immersion: Refractive index of the immersion medium

Returns

  • Apodization factor for energy conservation
source
MicroscopePSFs.calculate_axial_phaseFunction
calculate_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 medium
  • kz_coverslip: z-component of wave vector in coverslip
  • kz_immersion: z-component of wave vector in immersion medium
  • coverslip_thickness: Thickness of coverslip in mm (default: 0.17mm)

Returns

  • Defocus phase
source
MicroscopePSFs.calculate_dipole_field_componentsMethod
calculate_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 plane
  • sinθ: Sine of polar angle
  • cosθ: Cosine of polar angle
  • dipole: Dipole orientation vector
  • Tp: p-polarization transmission coefficient
  • Ts: s-polarization transmission coefficient

Returns

  • Tuple (Ex, Ey) of complex field components
source
MicroscopePSFs.calculate_interface_fresnelMethod
calculate_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 microns
  • n1: Refractive index of first medium
  • n2: Refractive index of second medium

Returns

  • Tuple of (Tp, Ts): p and s polarization transmission coefficients
source
MicroscopePSFs.calculate_wave_vectorsMethod
calculate_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 microns
  • n_medium: Refractive index of the sample medium
  • n_coverslip: Refractive index of the coverslip
  • n_immersion: Refractive index of the immersion medium

Returns

  • Tuple of (kzmedium, kzcoverslip, kz_immersion): z-components of wave vectors in each medium
source
MicroscopePSFs.fill_vector_pupils!Function
fill_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 fill
  • dipole: Dipole orientation vector
  • base_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)
source
MicroscopePSFs.get_pixel_indicesMethod
get_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 coordinates
  • emitter: Emitter with position information
  • support: 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
source
MicroscopePSFs.integrate_pixels!Method
integrate_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 stored
  • psf: Point spread function to integrate
  • camera: Camera geometry defining pixel edges
  • emitter: Emitter with position information
  • sampling: Subpixel sampling density (default: 2)
  • threaded: Whether to use multi-threading for integration (default: true)

Returns

  • The result array, now filled with integrated intensities
source
MicroscopePSFs.integrate_pixels!Method
integrate_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 stored
  • psf: Point spread function to integrate
  • pixel_edges_x, pixel_edges_y: Arrays defining pixel edge coordinates
  • emitter: Emitter with position information
  • 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

  • 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
source
MicroscopePSFs.integrate_pixels_amplitude!Method
integrate_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 stored
  • psf: Point spread function to integrate
  • camera: Camera geometry defining pixel edges
  • emitter: Emitter with position information
  • sampling: 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
source
MicroscopePSFs.integrate_pixels_amplitude!Method
integrate_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 stored
  • psf: Point spread function to integrate
  • pixel_edges_x, pixel_edges_y: Arrays defining pixel edge coordinates
  • emitter: Emitter with position information
  • sampling: 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
source
MicroscopePSFs.integrate_pixels_amplitude!Method
integrate_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 instance
  • camera: Camera geometry defining pixel edges
  • emitter: Emitter with position information
  • sampling: 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
source
MicroscopePSFs.normalize!Method
normalize!(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
source
MicroscopePSFs.update_pupils!Method
update_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
source

Zernike Module Internal API