Zernike.plotconfig — Constant
Zernike plot settings.
Fields / Options:
size::Tuple{Float64, Float64}: window size (DPI scaled resolution);fontsize::Float64: text size;colormap::Symbol: Default::oslo;theme::Makie.Attributes: Default:theme_black();focus_on_show::Bool: whether the window is focused on generation (default:true).
There are two methods which can be used to trigger a settings refresh: resize! & reset!.
The focus_on_show attribute of plotconfig only controls whether Zernike plots will set that option of the screen configuration, but if using Makie to create other types of plots this property needs to be set using activate! or set_theme!. The Zernike.reset!(plotconfig) function will appropriately reset these along with the rest of the plotconfig.
See also: zplot.
Zernike.Gradient — Type
Zernike.Gradient(Z::Polynomial)Returns the gradient ∇Z(ρ, θ) in a polar basis.
If called with a complex number x + iy this function returns the vector [∂Z/∂x, ∂Z/∂y] at the point (x, y) in Cartesian coordinates instead.
See also: grad, lap, derivatives, Laplacian.
Zernike.Laplacian — Type
Zernike.Laplacian(Z::Polynomial)Returns the Laplacian ΔZ(ρ, θ).
This function can be evaluated in Cartesian coordinates if passed a complex number instead (i.e. ∇²Z(x, y) ≡ ΔZ(xy::Complex)).
See also: grad, lap, derivatives, Gradient.
Zernike.Polynomial — Type
Zernike.Polynomial
Callable type: function Z(ρ, θ) bound to a given set of Zernike indices m and n.
The single argument method Zₘₙ(ρ) will radially evaluate the polynomial with angle zero.
Fields:
inds: named tuple containing the Zernike polynomial indices;N: normalization factor;R:RadialPolynomialcallable type: functionR(ρ);M:Harmoniccallable type: functionM(θ).
This type can be indexed (zero-based) to return a specific radial coefficient corresponding to the term with exponent i. Calling getindex without an explicit index will return the full vector of coefficients.
See also: Zernike.Wavefront.
Zernike.Standard — Type
Standard(v::Vector{Float64})Wraps a standard ANSI / ISO vector of Zernike polynomial single indices for inter-conversions, viz. Noll(s::Standard) & Fringe(s::Standard).
See also: standardize.
Zernike.Wavefront — Type
Zernike.Wavefront
Callable type: function ΔW(ρ, θ) bound to a given set of Zernike polynomial functions Zᵢ(ρ, θ) and their corresponding expansion coefficients aᵢ.
Specifically, ΔW(ρ, θ) = ∑aᵢZᵢ(ρ, θ)
The single argument method ΔW(ρ) will radially evaluate the polynomials with angle zero.
Fields:
recap: vector of named tuples containing the Zernike polynomial indices and the corresponding expansion coefficients;v: full vector of the unfiltered full-precision standardized expansion coefficients up ton_max;n_max: maximum radial degree fit to;fit_to: vector of(m, n)tuples specifying the polynomials used for the fit;a: vector of the Zernike expansion coefficients corresponding to each polynomial present;Z: the respective Zernike polynomial functions;precision: the precision with which thea&Zvalues were determined;ssr: the sum of the squared residuals from the fit.
The fit_to field is an empty vector if the default full range up to n_max (0:j_max) was used with no orders specified. Note that these orders could differ from the polynomials determined after the fit; they are simply what was passed to the fitting function and may refer to polynomials not present in the reconstruction if after filtering the corresponding coefficients are zero.
This type can be indexed (zero-based) to return a specific Zernike expansion coefficient corresponding to the Zernike polynomial of index j. Calling getindex without an explicit index will return the full vector of coefficients.
See also: Zernike.Polynomial.
Wavefront(G::Gradient{Polynomial})Return the Zernike coefficients of the gradient as a Wavefront.
See also: grad.
Wavefront(L::Laplacian)Return the Zernike coefficients of the Laplacian as a Wavefront.
See also: lap.
Wavefront{RadialPolynomial}(m::Int, a::Vector)Construct a radial series in a Wavefront from an azimuthal index and a coefficient vector associated with a sequence of radial polynomials.
Wavefront{RadialPolynomial}(m::Int, n::Vector{Int}, a::Vector)Construct a Wavefront radial polynomial series from an azimuthal index, a vector of radial indices, and a coefficient vector associated with them; the indices are subject to the usual restrictions for Zernike polynomials.
See also *.
Wavefront(aberr::Aberration)Convert a set of Seidel aberrations to Zernike polynomials.
Base.:* — Function
*(w1::Wavefront{T}, w2::Wavefront{T}; [threshold = 5.0E-324]) where T <: RadialPolynomial
w1 * w2Compute the radial product expansion of two Zernike radial polynomial sequences.
Returns the coefficients of the new series embedded in a Wavefront. Based on a paper by Cadiot et al. (2024) using Matrix Multiplication Transforms and their inverses to perform quadrature.
The coefficients will be filtered according to threshold, with any absolute values less than it being set to zero.
See also: Wavefront.
Zernike.PSF — Function
PSF(ΔW::Wavefront)Return the incoherent Point Spread Function as a real valued matrix for the input Wavefront.
The output is normalized by the diffraction limited peak intensity i.e. the maximum is the Strehl ratio.
The keyword argument s (default 50.0) determines the scale, i.e. the Fourier space step size is 0.25 / s.
See also: See also: OTF, MTF, psf_plot, psf_plot!, mtf_plot, mtf_plot!.
Zernike.W — Function
W(ρ, θ, OPD, n_max)Return the Wavefront function ΔW(ρ, θ) corresponding to an n_max fit.
W(∂x::Vector{Float64}, ∂y::Vector{Float64}; [normalized = true])Compute the original wavefront error coefficients given the Zernike expansion coefficient vectors of the partial derivatives of the wavefront error.
normalized refers to whether the derivatives are associated with unnormalized Zernike polynomials.
W(∂x::Wavefront, ∂y::Wavefront)Equivalent to the above, but for input wavefronts. Useful if the derivative data was fit into a Zernike basis using W. Returns a Wavefront.
Zernike.format_strings — Function
format_strings(Z::AbstractPolynomial)Return a 3-tuple with the index formatted LaTeX variable name, the full LaTeX string equation, and the Unicode string representation of the polynomial. Calling the polynomial with the String type is a shortcut which yields just the symbolic Unicode string.
julia> Z(-8, 8)(String)
"√(18)ρ⁸sin(8θ)"format_strings(ΔW::Wavefront)Return a (possibly truncated) symbolic representation of the wavefront error in a Zernike basis as a LaTeXString.
See also: print_strings.
Zernike.fringe_to_j — Function
fringe_to_j(fringe::Int)Convert Fringe indices to ANSI standard indices.
Only indices 1:37 are valid.
See also: j_to_fringe, noll_to_j, j_to_noll, standardize, get_j, get_mn.
Zernike.get_j — Function
get_j(m::Int, n::Int)Return the single mode-ordering index j corresponding to azimuthal & radial indices (m, n).
get_j(n_max::Int)Return the single mode-ordering index j corresponding to the maximum radial index n_max; equivalent to get_j(n_max, n_max).
See also: get_mn, noll_to_j, j_to_noll, fringe_to_j, j_to_fringe, standardize.
Zernike.get_m — Function
get_m(j::Int)Return the azimuthal index m given the single mode-ordering index j.
See also: get_n, get_mn, get_j, noll_to_j, j_to_noll, fringe_to_j, j_to_fringe, standardize.
Zernike.get_mn — Function
get_mn(j::Int)Return the azimuthal & radial indices (m, n) given the single mode-ordering index j.
See also: get_j, noll_to_j, j_to_noll, fringe_to_j, j_to_fringe, standardize.
Zernike.get_n — Function
get_n(j::Int)Return the radial index n given the single mode-ordering index j.
See also: get_m, get_mn, get_j, noll_to_j, j_to_noll, fringe_to_j, j_to_fringe, standardize.
Zernike.grad — Function
∇(Z)
grad(Z::Polynomial)Return the gradient of the Zernike polynomial expressed as a Wavefront in Zernike polynomial expansion coefficients.
Equivalent to:
Wavefront(g::Gradient{Polynomial})
Wavefront(::Type{<:Gradient}, m, n; [normalize = true])
Wavefront(::Type{<:Gradient}, j; [normalize = true])with the last two methods allowing unnormalized input. Normalized means the polynomial is expressed as Z(ρ, θ) = N * R(ρ) * M(θ) with N being the normalization prefactor required so that π-normalized integration with respect to the areal measure of Z² over the unit disk yields unity .
grad(m, n, ::Type{Matrix{Complex}})Return a complex matrix encoding the partial derivatives of a complex unnormalized Zernike polynomial at indices m & n.
The first two columns of the matrix refer to the x partial derivatives of Z(|m|, n) & Z(-|m|, n), respectively. Likewise, the last two columns refer to the y partial derivatives of Z(|m|, n) & Z(-|m|, n).
grad(m, n, ::Type{Vector{Complex}})Convenience method which returns the relevant sign-dependent complex gradient from the above mentioned grad(m, n, Matrix{Complex}) call.
grad(m, n, ::Type{Vector{Real}}; [normalize = true])Return the gradient of the Zernike polynomial as a tuple of real-valued Zernike coefficient vectors.
∇(W)
grad(W::Wavefront)Computes the gradient of the wavefront.
See also: Zernike.lap, Zernike.Gradient, Zernike.Laplacian, derivatives.
Zernike.j_to_fringe — Function
j_to_fringe(j::Int)Convert ANSI standard indices to Fringe indices.
Call fringe_to_j.(1:37) to return valid indices.
See also: fringe_to_j, noll_to_j, j_to_noll, standardize, get_j, get_mn.
Zernike.j_to_noll — Function
j_to_noll(j::Int)Convert ANSI standard indices to Noll indices.
See also: noll_to_j, fringe_to_j, j_to_fringe, get_j, get_mn, standardize.
Zernike.lap — Function
lap(Z::Polynomial)Return the Laplacian of the Zernike polynomial expressed as a Wavefront in Zernike polynomial expansion coefficients.
Equivalent to:
Wavefront(l::Laplacian)
Wavefront(::Type{<:Laplacian}, m, n; [normalize = true])
Wavefront(::Type{<:Laplacian}, j; [normalize = true])with the last two methods allowing unnormalized input. Normalized means the polynomial is expressed as Z(ρ, θ) = N * R(ρ) * M(θ) with N being the normalization prefactor required so that π-normalized integration with respect to the areal measure of Z² over the unit disk yields unity .
lap(m, n; [normalize = true])Return the Laplacian of the Zernike polynomial as a real-valued coefficient vector.
See also: grad, Gradient, Laplacian, derivatives.
Zernike.map_phase — Function
map_phase(ρ, θ, OPD)Reverse dimensional coordinate transform with respect to the main wavefront error method. Returns the OPD as a matrix along with the corresponding unique coordinate vectors. Assumes uniform sampling.
See also: wavefront.
Zernike.metrics — Function
metrics(ΔW::Wavefront)Compute wavefront error metrics. Returns a named 3-tuple with the peak-to-valley error, RMS wavefront error, and Strehl ratio.
Zernike.mnv — Function
mnv(v)Zernike tool which takes a vector v (e.g. a vector of single-index j integers or a real or complex floating point coefficient vector) and produces a (length(v) × 3) Matrix{Number} with columns corresponding to the indices m, n, & the vector v.
Zernike.mtf_plot — Function
mtf_plot(mtf::Matrix)Plot the Modulation Transfer Function as a three-dimensional surface plot.
mtf_plot(mtf::Matrix, x_or_y::Symbol)Plot the MTF along the sagittal or tangential direction in the pupil. The input symbol must be either :x or :y corresponding to the relevant meridian.
See also: mtf_plot!, MTF, OTF, zplot, plotconfig, PSF, psf_plot, psf_plot!.
Zernike.mtf_plot! — Function
mtf_plot!(mtf::Matrix, x_or_y::Symbol)Overlay the MTF on the existing plot. The input symbol must be either :x or :y corresponding to the relevant meridian.
See also: mtf_plot, MTF, OTF, zplot, plotconfig, PSF, psf_plot, psf_plot!.
Zernike.noll_to_j — Function
noll_to_j(noll::Int)Convert Noll indices to ANSI standard indices.
See also: j_to_noll, fringe_to_j, j_to_fringe, get_j, get_mn, standardize.
Zernike.print_strings — Function
print_strings([io::IO], j::AbstractVector; [limit = true])Print the Unicode string representations of select Zernike polynomials; j here can be a range / vector of Zernike single indices; the output stream defaults to stdout and the vertical output is limited by your displaysize unless you pass the keyword argument limit = false.
print_strings(j_max::Int)Print the Unicode string representations of the first j + 1 Zernike polynomials from single index 0 to j_max.
See also: format_strings.
Zernike.psf_plot — Function
psf_plot(psf::Matrix)Plot the Point Spread Function as a three-dimensional surface plot.
psf_plot(psf::Matrix, x_or_y::Symbol)Plot the PSF along the sagittal or tangential direction. The input symbol must be either :x or :y corresponding to the relevant meridian.
See also: psf_plot!, PSF, mtf_plot, mtf_plot!, MTF, OTF, zplot, plotconfig.
Zernike.psf_plot! — Function
psf_plot!(psf::Matrix, x_or_y::Symbol)Overlay the PSF on the existing plot. The input symbol must be either :x or :y corresponding to the relevant meridian.
See also: psf_plot, PSF, mtf_plot, mtf_plot!, MTF, OTF, zplot, plotconfig.
Zernike.radial_coefficients — Function
radial_coefficients(m::Int, n::Int, T::Type{<:Number} = Float64)Compute Zernike radial polynomial coefficients as type T using an algorithm based on Honarvar & Paramesran's recursive relation suitable for high orders.
The coefficients in the vector correspond to terms with powers in ascending order for a Zernike polynomial with indices m & n subject to the following requirements:
n ≥ 0
|m| ≤ n
n ≡ m (mod 2).See also: wavefront_coefficients, transform_coefficients.
Zernike.reduce_wave — Function
reduce_wave(W::Wavefront, precision::Int)Reduces Wavefront precision.
Zernike.reset! — Function
reset!(plotconfig::PlotConfig; [only_theme = false])Reset all of the Zernike.plotconfig settings to their defaults.
This will also reset the part of the GLMakie global / universal theme used by Zernike, particularly the window title and focus_on_show properties, until another Zernike plot is crafted. If you only want to reset the theme and keep the plotconfig intact then call with the keyword argument only_theme = true.
Zernike.resize! — Function
resize!(plotconfig::PlotConfig)Reset only the size and fontsize settings for Zernike.plotconfig. This is useful if your primary monitor changes or you want to return to the automatically determined values.
Zernike.scale — Function
scale(v, ε; precision, finesse)Scale the pupil over a wavefront using an algorithm based on Janssen & Dirksen's formula and plot the result.
v is the set of Zernike wavefront error expansion coefficients and ε is the scaling factor.
Zernike.sieve — Function
sieve(v::Vector{Float64}, threshold::Float64)Zero out any elements lower than the threshold.
sieve(a::Vector)Return Zernike non-sequential indices and the corresponding wavefront expansion coefficients for any non-zero coefficients in the input vector.
Zernike.standardize — Function
standardize(noll::Noll)
standardize(fringe::Fringe)Format a Noll or Fringe specified Zernike expansion coefficient vector according to the ANSI standard.
Floating-point coefficient vectors need to be wrapped in the index types (e.g. standardize(Fringe(v))).
The Fringe method expects unnormalized coefficients; the input coefficients will be re-ordered and normalized in line with the orthonormal standard. As Fringe is a 37 polynomial subset of the full set of Zernike polynomials any coefficients in the standard order missing a counterpart in the input vector will be set to zero.
See also: noll_to_j, j_to_noll, fringe_to_j, j_to_fringe, get_j, get_mn.
standardize(v_sub::FloatVec, [j::AbstractVector{Int}])
standardize(v_sub::Vector, orders::Vector{Tuple{Int, Int}})
standardize(W::Wavefront)Pad a subset Zernike expansion coefficient vector to the full standard length up to n_max (1:j_max+1).
The tuples in orders must be of the form (m, n) associated with the respective coefficients at each index in v_sub.
j is a vector of single-mode ordering indices associated with the coefficients; if this is not supplied the coefficients will be assumed to be in order (0:j).
The Wavefront method pads the W.a coefficient vector.
Zernike.transform — Function
transform(v::Vector{T}, ε::T, [δ::Complex{T}], [ϕ::T], [ω::Tuple{T,T}]) where T <: Float64Compute a new set of Zernike wavefront error expansion coefficients under a given set of transformation factors and plot the result.
Available transformations are scaling, translation, & rotation for circular and elliptical exit pupils. These are essentially coordinate transformations in the pupil plane over the wavefront map.
Main arguments
v: vector of full Zernike expansion coefficients ordered in accordance with the ANSI / OSA single index standard. This is thevvector returned bywavefront(ρ, θ, OPD, n_max);ε: scaling factor{0 ≤ ε ≤ 1};δ: translational complex coordinates (displacement of the pupil center in the complex plane);ϕ: rotation of the pupil in radians(mod 2π), defined positive counter-clockwise from the horizontal x-axis;ω: elliptical pupil transform parameters; 2-tuple whereω[1]is the ratio of the semi-minor axis length to the semi-major axis length of the ellipse andω[2]is the angle defined positive counter-clockwise from the horizontal coordinate axis of the exit pupil to the minor axis of the ellipse.
The order the transformations are applied is:
scaling –> translation –> rotation –> elliptical transform.
See also: Y, zernike, wavefront, transform_coefficients.
Keyword argument options:
transform(v, ε, [δ], [ϕ], [ω]; [precision = 3], [finesse::Int])precision: number of digits to use after the decimal point in computing the expansion coefficients. Results will be rounded according to this precision and any polynomials with zero-valued coefficients will be ignored when pulling in the Zernike functions while constructing the composite wavefront error; this means lower precision values yield faster results.finesse: determines the size of the plot matrix; sets both dimensions to the given value.
Extended help
ε = r₂/r₁ where r₂ is the new smaller radius, r₁ the original
In particular the radial variable corresponding to the rescaled exit pupil is normalized such that:ρ = r/r₂; {0 ≤ ρ ≤ 1}r: radial pupil position, r₂: max. radiusΔW₂(ρ₂, θ) = ΔW₁(ερ₂, θ)
For translation the shift must be within the bounds of the scaling applied such that:0.0 ≤ ε + |δ| ≤ 1.0.
For elliptical pupils (usually the result of measuring the wavefront off-axis), the semi-major axis is defined such that it equals the radius of the circle and so ω[1] is the fraction of the circular pupil covered by the semi-minor axis (this is approximated well by a cosine projection factor for angles up to 40 degrees); ω[2] is then the direction of the stretching applied under transformation in converting the ellipse to a circle before fitting the expansion coefficients.
The transformed expansion coefficients are computed using a fast and accurate algorithm suitable for high orders; it is based on a formulation presented by Lundström & Unsbo (2007) doi:10.1364/JOSAA.24.000569.
Zernike.transform_coefficients — Function
transform_coefficients(v, ε, δ, ϕ, ω)Directly compute Zernike wavefront error expansion coefficients under pupil transformations. The argument types are the same as in transform.
Returns a 2-tuple with the new coefficient vector and order n_max.
See also: transform, radial_coefficients, wavefront_coefficients.
Zernike.wavefront — Function
wavefront(ρ, θ, OPD, n_max)Fit wavefront errors up to order n_max.
Estimates wavefront error by expressing optical aberrations as a linear combination of weighted Zernike polynomials using a linear least squares method. The accuracy of this type of wavefront reconstruction represented as an expanded series depends upon a sufficiently sampled phase field and a suitable choice of the fitting order n_max.
Main arguments
ρ, θ, and OPD must be floating-point vectors of equal length; at each specific index the values are elements of an ordered triple over the exit pupil.
ρ: normalized radial exit pupil position variable{0 ≤ ρ ≤ 1};θ: angular exit pupil variable in radians(mod 2π), defined positive counter-clockwise from the horizontal x-axis;OPD: measured optical path difference in waves;n_max: maximum radial degree to fit to.
Return values
Returns six values contained within a WavefrontOutput type, with fields:
recap: vector of named tuples containing the Zernike polynomial indices and the corresponding expansion coefficients rounded according toprecision;v: full vector of Zernike wavefront error expansion coefficients;ssr: the sum of the squared residuals from the fit;metrics: named 3-tuple with the peak-to-valley error, RMS wavefront error, and Strehl ratio;W: theWavefrontfunctionΔW(ρ, θ);fig: theMakieFigureAxisPlot.
See also: W, zernike, transform.
wavefront(ρ, θ, OPD, orders::Vector{Tuple{Int, Int}})Fit wavefront errors to specific Zernike polynomials specified in orders containing Zernike (m, n) tuples.
wavefront(OPD, fit_to; options...)Fitting method accepting a floating-point matrix of phase data uniformly produced in a polar coordinate system over the pupil.
The matrix is expected to be a polar grid of regularly spaced periodic samples with the first element referring to the value at the origin and the end points including the boundary of the pupil (i.e. ρ, θ = 0.0:step:1.0, 0.0:step:2π). The first axis of the matrix (the rows) must correspond to the angular variable θ while the second axis (the columns) must correspond to the radial variable ρ.
fit_to can be either n_max::Int or orders::Vector{Tuple{Int, Int}}.
wavefront(ρ::Vector, θ::Vector, OPD::Matrix, fit_to; options...)Fitting method accepting coordinate vectors and a floating-point matrix of corresponding phase data produced in a polar coordinate system over the pupil under the aforementioned dimensional ordering assumption. This method does not assume equally spaced samples.
wavefront(x, y, OPD; fit_to, options...)Fitting method accepting normalized Cartesian coordinate data.
If OPD is a matrix the shape of the axes is assumed to be x-by-y.
Keyword argument options:
wavefront(ρ, θ, OPD, n_max; [precision = 3], [finesse::Int])precision: number of digits to use after the decimal point in computing the expansion coefficients. Results will be rounded according to this precision and any polynomials with zero-valued coefficients will be ignored when pulling in the Zernike functions while constructing the composite wavefront error; this means lower precision values yield faster results.finesse: determines the size of the plot matrix; sets both dimensions to the given value.
Zernike.wavefront_coefficients — Function
wavefront_coefficients(ρ, θ, OPD, n_max)Returns a 2-tuple with the full vector of Zernike expansion coefficients obtained through the least squares fit and the corresponding sum of the squared residuals.
See also: wavefront, radial_coefficients, transform_coefficients.
Zernike.zernike — Function
zernike(m, n)
zernike(j)Plot a Zernike polynomial of azimuthal order m and radial degree n.
The single index j begins at zero and follows the ANSI Z80.28-2004 / ISO 24157:2008 / Optica (OSA) standard.
Returns a Zernike.Output type which contains (among other things):
Z: thePolynomialfunctionZ(ρ, θ);fig: theMakieFigureAxisPlot;coeffs: vector of radial polynomial coefficients;latex:LaTeXstring of the Zernike polynomial;unicode:Unicodestring of the Zernike polynomial.
The coefficients belong to terms with exponent n - 2(i - 1) where i is the vector's index.
The radial polynomial coefficients are computed using a fast and accurate algorithm suitable for high orders; it is based on a recursive relation presented by Honarvar & Paramesran (2013) doi:10.1364/OL.38.002487.
See also: Z, wavefront, transform, radial_coefficients.
Keyword argument options:
zernike(m, n; [finesse::Int = 1024])finesse: determines the size of the plot matrix; sets both dimensions to the given value.
Zernike.zplot — Function
zplot(φ; kwargs...)
zplot(ρ, θ, φ; kwargs...)Plot Zernike phase function types (Polynomials, Wavefronts, Derivatives, arithmetic types, etc.) as well as quantized phase arrays; for the latter the arguments must be a collection of discretized samples where the polar variable objects refer to either ranges, vectors, or 1-dimensional matrices & the like, and the phase structure φ is either an array corresponding to these samples or a callable type in which case the matrix will be constructed for you.
Keyword arguments:
size::Tuple{Float64, Float64}: window size (DPI scaled resolution);fontsize::Float64: text size;colormap::Symbol: Default::oslo;theme::Makie.Attributes: Default:theme_black();focus_on_show::Bool: whether the window is focused on generation (default:true);window_title::String: window title;plot_title::Union{String, LaTeXString}: plot title;m::Int: azimuthal order (used to determine matrix size);n::Int: radial order (used to determine matrix size);finesse::Int: determines the matrix size;high_order::Bool: whether to apply a logarithmic transform (default:false).
Any keyword arguments supported by Makie's surface are also supported.
Plots can be updated on demand by passing an Observable and changing its value.
For example:
w = Observable(Wavefront([0.0, -1.0, 1.0]))
zplot(w)
# update
w[] = Wavefront([0.0, 1.0, 1.0])As a convenient shortcut any type of phase object can be plotted by simply calling it with no arguments, e.g. as w(); similarly, calling it as w(Screen) will plot it in a new window; note the first method depends on display automatically being called, while the second will explicitly call it.
See also: plotconfig.
Zernike.∇ — Function
∇(Z)
grad(Z::Polynomial)Return the gradient of the Zernike polynomial expressed as a Wavefront in Zernike polynomial expansion coefficients.
Equivalent to:
Wavefront(g::Gradient{Polynomial})
Wavefront(::Type{<:Gradient}, m, n; [normalize = true])
Wavefront(::Type{<:Gradient}, j; [normalize = true])with the last two methods allowing unnormalized input. Normalized means the polynomial is expressed as Z(ρ, θ) = N * R(ρ) * M(θ) with N being the normalization prefactor required so that π-normalized integration with respect to the areal measure of Z² over the unit disk yields unity .
grad(m, n, ::Type{Matrix{Complex}})Return a complex matrix encoding the partial derivatives of a complex unnormalized Zernike polynomial at indices m & n.
The first two columns of the matrix refer to the x partial derivatives of Z(|m|, n) & Z(-|m|, n), respectively. Likewise, the last two columns refer to the y partial derivatives of Z(|m|, n) & Z(-|m|, n).
grad(m, n, ::Type{Vector{Complex}})Convenience method which returns the relevant sign-dependent complex gradient from the above mentioned grad(m, n, Matrix{Complex}) call.
grad(m, n, ::Type{Vector{Real}}; [normalize = true])Return the gradient of the Zernike polynomial as a tuple of real-valued Zernike coefficient vectors.
∇(W)
grad(W::Wavefront)Computes the gradient of the wavefront.
See also: Zernike.lap, Zernike.Gradient, Zernike.Laplacian, derivatives.