Analysis

The batoid analysis module contains functions for evaluating wavefronts, PSFs, and Zernike polynomial decompositions.

batoid.analysis.dkdu(optic, theta_x, theta_y, wavelength, nrad=6, naz=36, projection='postel')[source]

Calculate derivative of outgoing ray k-vector with respect to incoming ray pupil coordinate.

Parameters:
  • optic (batoid.Optic) – Optical system

  • theta_x, theta_y (float) – Field angle in radians

  • wavelength (float) – Wavelength in meters

  • nrad (int, optional) – Number of ray radii to use. (see RayVector.asPolar())

  • naz (int, optional) – Approximate number of azimuthal angles in outermost ring. (see RayVector.asPolar())

  • projection ({‘postel’, ‘zemax’, ‘gnomonic’, ‘stereographic’, ‘lambert’, ‘orthographic’}) – Projection used to convert field angle to direction cosines.

Returns:

dkdu ((2, 2), ndarray) – Jacobian transformation matrix for converting between (kx, ky) of rays impacting the focal plane and initial pupil plane coordinate.

batoid.analysis.drdth(optic, theta_x, theta_y, wavelength, nrad=6, naz=36, projection='postel')[source]

Calculate derivative of focal plane coord with respect to field angle.

Parameters:
  • optic (batoid.Optic) – Optical system

  • theta_x, theta_y (float) – Field angle in radians

  • wavelength (float) – Wavelength in meters

  • nrad (int, optional) – Number of ray radii to use. (see RayVector.asPolar())

  • naz (int, optional) – Approximate number of azimuthal angles in outermost ring. (see RayVector.asPolar())

  • projection ({‘postel’, ‘zemax’, ‘gnomonic’, ‘stereographic’, ‘lambert’, ‘orthographic’}) – Projection used to convert field angle to direction cosines.

Returns:

drdth ((2, 2), ndarray) – Jacobian transformation matrix for converting between (theta_x, theta_y) and (x, y) on the focal plane.

Notes

This is the Jacobian of pixels -> tangent plane, (and importantly, not pixels -> ra/dec). It should be close to the inverse plate scale though, especially near the center of the tangent plane projection.

batoid.analysis.dthdr(optic, theta_x, theta_y, wavelength, nrad=6, naz=36, projection='postel')[source]

Calculate derivative of field angle with respect to focal plane coordinate.

Parameters:
  • optic (batoid.Optic) – Optical system

  • theta_x, theta_y (float) – Field angle in radians

  • wavelength (float) – Wavelength in meters

  • nrad (int, optional) – Number of ray radii to use. (see RayVector.asPolar())

  • naz (int, optional) – Approximate number of azimuthal angles in outermost ring. (see RayVector.asPolar())

  • projection ({‘postel’, ‘zemax’, ‘gnomonic’, ‘stereographic’, ‘lambert’, ‘orthographic’}) – Projection used to convert field angle to direction cosines.

Returns:

dthdr ((2, 2), ndarray) – Jacobian transformation matrix for converting between (x, y) on the focal plane and field angle (theta_x, theta_y).

Notes

This is the Jacobian of tangent plane -> pixels, (and importantly, not ra/dec -> pixels). It should be close to the plate scale though, especially near the center of the tangent plane projection.

batoid.analysis.huygensPSF(optic, theta_x, theta_y, wavelength, projection='postel', nx=None, dx=None, dy=None, nxOut=None, reference='chief')[source]

Compute a PSF via the Huygens construction.

Parameters:
  • optic (batoid.Optic) – Optical system

  • theta_x, theta_y (float) – Field angle in radians

  • wavelength (float) – Wavelength in meters

  • projection ({‘postel’, ‘zemax’, ‘gnomonic’, ‘stereographic’, ‘lambert’, ‘orthographic’}) – Projection used to convert field angle to direction cosines.

  • nx (int, optional) – Size of ray grid to use.

  • dx, dy (float, optional) – Lattice scales to use for PSF evaluation locations. Default, use fftPSF lattice.

  • nxOut (int, optional) – Size of the output lattice. Default is to use nx.

  • reference ({‘chief’, ‘mean’}) – If ‘chief’, then center the output lattice where the chief ray intersects the focal plane. If ‘mean’, then center at the mean non-vignetted ray intersections.

Returns:

psf (batoid.Lattice) – The PSF.

Notes

The Huygens construction is to evaluate the PSF as

\[I(x) \propto \sum_u \exp(i \phi(u)) \exp(i k(u) \cdot r)\]

The \(u\) are assumed to uniformly sample the entrance pupil, but not include any rays that get vignetted before they reach the focal plane. The \(\phi\) s are the phases of the exit rays evaluated at a single arbitrary time. The \(k(u)\) indicates the conversion of the uniform entrance pupil samples into nearly (though not exactly) uniform samples in k-space of the output rays.

The output locations where the PSF is evaluated are governed by dx, dy, and nx. If dx and dy are None, then the same lattice as in fftPSF will be used. If dx and dy are scalars, then a lattice with primitive vectors [dx, 0] and [0, dy] will be used. If dx and dy are 2-vectors, then those will be the primitive vectors of the output lattice.

batoid.analysis.wavefront(optic, theta_x, theta_y, wavelength, projection='postel', nx=32, sphereRadius=None, reference='chief')[source]

Compute wavefront.

Parameters:
  • optic (batoid.Optic) – Optical system

  • theta_x, theta_y (float) – Field angle in radians

  • wavelength (float) – Wavelength in meters

  • projection ({‘postel’, ‘zemax’, ‘gnomonic’, ‘stereographic’, ‘lambert’, ‘orthographic’}) – Projection used to convert field angle to direction cosines.

  • nx (int, optional) – Size of ray grid to use.

  • sphereRadius (float, optional) – The radius of the reference sphere. Nominally this should be set to the distance to the exit pupil, though the calculation is usually not very sensitive to this. Many of the telescopes that come with batoid have values for this set in their yaml files, which will be used if this is None.

  • reference ({‘chief’, ‘mean’}) – If ‘chief’, then center the output lattice where the chief ray intersects the focal plane. If ‘mean’, then center at the mean non-vignetted ray intersections.

Returns:

wavefront (batoid.Lattice) – A batoid.Lattice object containing the wavefront values in waves and the primitive lattice vectors of the entrance pupil grid in meters.

batoid.analysis.fftPSF(optic, theta_x, theta_y, wavelength, projection='postel', nx=32, pad_factor=2, sphereRadius=None, reference='chief')[source]

Compute PSF using FFT.

Parameters:
  • optic (batoid.Optic) – Optical system

  • theta_x, theta_y (float) – Field angle in radians

  • wavelength (float) – Wavelength in meters

  • projection ({‘postel’, ‘zemax’, ‘gnomonic’, ‘stereographic’, ‘lambert’, ‘orthographic’}) – Projection used to convert field angle to direction cosines.

  • nx (int, optional) – Size of ray grid to use.

  • pad_factor (int, optional) – Factor by which to pad pupil array. Default: 2

  • sphereRadius (float, optional) – The radius of the reference sphere. Nominally this should be set to the distance to the exit pupil, though the calculation is usually not very sensitive to this. Many of the telescopes that come with batoid have values for this set in their yaml files, which will be used if this is None.

  • reference ({‘chief’, ‘mean’}) – If ‘chief’, then center the output lattice where the chief ray intersects the focal plane. If ‘mean’, then center at the mean non-vignetted ray intersection.

Returns:

psf (batoid.Lattice) – A batoid.Lattice object containing the relative PSF values and the primitive lattice vectors of the focal plane grid.

batoid.analysis.zernike(optic, theta_x, theta_y, wavelength, projection='postel', nx=32, sphereRadius=None, reference='chief', jmax=22, eps=0.0)[source]

Compute Zernike polynomial decomposition of the wavefront.

This calculation propagates a square grid of rays to the exit pupil reference sphere where the wavefront is computed. The optical path differences of non-vignetted rays are then fit to Zernike polynomial coefficients numerically.

Parameters:
  • optic (batoid.Optic) – Optical system

  • theta_x, theta_y (float) – Field angle in radians

  • wavelength (float) – Wavelength in meters

  • projection ({‘postel’, ‘zemax’, ‘gnomonic’, ‘stereographic’, ‘lambert’, ‘orthographic’}) – Projection used to convert field angle to direction cosines.

  • nx (int, optional) – Size of ray grid to use.

  • sphereRadius (float, optional) – The radius of the reference sphere. Nominally this should be set to the distance to the exit pupil, though the calculation is usually not very sensitive to this. Many of the telescopes that come with batoid have values for this set in their yaml files, which will be used if this is None.

  • reference ({‘chief’, ‘mean’}) – If ‘chief’, then center the output lattice where the chief ray intersects the focal plane. If ‘mean’, then center at the mean non-vignetted ray intersection.

  • jmax (int, optional) – Number of coefficients to compute. Default: 12.

  • eps (float, optional) – Use annular Zernike polynomials with this fractional inner radius. Default: 0.0.

Returns:

zernikes (array) – Zernike polynomial coefficients in waves.

Notes

Zernike coefficients are indexed following the Noll convention. Additionally, since python lists start at 0, but the Noll convention starts at 1, the 0-th index of the returned array is meaningless. I.e., zernikes[1] is piston, zernikes[4] is defocus, and so on…

Also, since Zernike polynomials are orthogonal over a circle or annulus, but vignetting may make the actual wavefront region of support something different, the values of fit coefficients can depend on the total number of coefficients being fit. For example, the j=4 (defocus) coefficient may depend on whether jmax=11 and jmax=21 (or some other value). See the zernikeGQ function for an alternative algorithm that is independent of jmax.

batoid.analysis.zernikeGQ(optic, theta_x, theta_y, wavelength, projection='postel', rings=6, spokes=None, sphereRadius=None, reference='chief', jmax=22, eps=0.0)[source]

Compute Zernike polynomial decomposition of the wavefront.

This calculation uses Gaussian Quadrature points and weights to compute the Zernike decomposition of the wavefront. The wavefront values at the GQ points are determined by tracing from the entrance pupil to the exit pupil reference sphere, and evaluating the optical path differences on this sphere.

Parameters:
  • optic (batoid.Optic) – Optical system

  • theta_x, theta_y (float) – Field angle in radians

  • wavelength (float) – Wavelength in meters

  • projection ({‘postel’, ‘zemax’, ‘gnomonic’, ‘stereographic’, ‘lambert’, ‘orthographic’}) – Projection used to convert field angle to direction cosines.

  • rings (int, optional) – Number of Gaussian quadrature rings to use. Default: 6.

  • spokes (int, optional) – Number of Gaussian quadrature spokes to use. Default: 2*rings + 1

  • sphereRadius (float, optional) – The radius of the reference sphere. Nominally this should be set to the distance to the exit pupil, though the calculation is usually not very sensitive to this. Many of the telescopes that come with batoid have values for this set in their yaml files, which will be used if this is None.

  • reference ({‘chief’, ‘mean’}) – If ‘chief’, then center the output lattice where the chief ray intersects the focal plane. If ‘mean’, then center at the mean non-vignetted ray intersection.

  • jmax (int, optional) – Number of coefficients to compute. Default: 22.

  • eps (float, optional) – Use annular Zernike polynomials with this fractional inner radius. Default: 0.0.

Returns:

zernikes (array) – Zernike polynomial coefficients in waves.

Notes

Zernike coefficients are indexed following the Noll convention. Additionally, since python lists start at 0, but the Noll convention starts at 1, the 0-th index of the returned array is meaningless. I.e., zernikes[1] is piston, zernikes[4] is defocus, and so on…

This algorithm takes advantage of the fact that Zernike polynomials are orthogonal on a circle or annulus to compute the value of individual coefficients via an integral:

\[a_j \propto \int Z_j(x, y) W(x, y) dx dy\]

This integral is approximated using Gaussian quadrature. Since the above integral depends on orthogonality, the wavefront must be decomposed over a circular or annular region of support. As such, this algorithm is unaffected by vignetting. It is required that no rays fail to be traced, even in vignetted regions.

batoid.analysis.zernikeTA(optic, theta_x, theta_y, wavelength, projection='postel', nrad=6, naz=36, reference='chief', jmax=22, eps=0.0, focal_length=None)[source]

Compute Zernikes in such a way that transverse ray aberrations are proportional to the gradient of the wavefront.

Parameters:
  • optic (batoid.Optic) – Optical system

  • theta_x, theta_y (float) – Field angle in radians

  • wavelength (float) – Wavelength in meters

  • projection ({‘postel’, ‘zemax’, ‘gnomonic’, ‘stereographic’, ‘lambert’, ‘orthographic’}) – Projection used to convert field angle to direction cosines.

  • nrad (int, optional) – Number of ray radii to use. (see RayVector.asPolar())

  • naz (int, optional) – Approximate number of azimuthal angles in outermost ring. (see RayVector.asPolar())

  • reference ({‘chief’, ‘mean’}) – If ‘chief’, then center the output lattice where the chief ray intersects the focal plane. If ‘mean’, then center at the mean non-vignetted ray intersection.

  • jmax (int, optional) – Number of coefficients to compute. Default: 22.

  • eps (float, optional) – Use annular Zernike polynomials with this fractional inner radius. Default: 0.0.

  • focal_length (float, optional) – Value to use for focal length. If omitted, then compute focal length from dr/dth.

Returns:

zernikes (array) – Zernike polynomial coefficients in waves.

Notes

Zernike coefficients are indexed following the Noll convention. Additionally, since python lists start at 0, but the Noll convention starts at 1, the 0-th index of the returned array is meaningless. I.e., zernikes[1] is piston, zernikes[4] is defocus, and so on…

batoid.analysis.doubleZernike(optic, field, wavelength, rings=6, spokes=None, kmax=22, **kwargs)[source]

Compute double Zernike polynomial decomposition of the wavefront.

The double Zernike decomposition describes both the focal and pupil variation of the wavefront. More specifically:

\[W(u, \theta) = \sum_{jk} a_{jk} Z_j(u), Z_k(\theta)\]

where \(u\) indicates the pupil coordinate and \(\theta\) indicates the field angle coordinate.

Parameters:
  • optic (batoid.Optic) – Optical system

  • field (float) – Outer field angle radius in radians.

  • wavelength (float) – Wavelength in meters

  • rings (int, optional) – Number of Gaussian quadrature rings to use. Default: 6.

  • spokes (int, optional) – Number of Gaussian quadrature spokes to use. Default: 2*rings + 1

  • kmax (int, optional) – Number of focal coefficients to compute. Default: 22.

  • **kwargs (dict) – Keyword arguments to pass to zernikeGQ.

Returns:

dzs (array) – Double Zernike polynomial coefficients in waves.

batoid.analysis.exitPupilPos(optic, wavelength, smallAngle=4.84813681109536e-06, **kwargs)[source]

Compute position of the exit pupil.

Traces a collection of small angle chief rays into object space, and then finds the mean of their closest approaches to the optic axis. Possibly only accurate for axially symmetric optical systems.

Parameters:
  • optic (batoid.Optic) – Optical system

  • wavelength (float) – Wavelength in meters

  • smallAngle (float, optional) – Angle in radians from which to search for the exit pupil position.

Returns:

Location of exit pupil in global coordinates (in meters).

batoid.analysis.focalLength(*args, **kwargs)[source]

Calculate focal length.

Parameters:
  • optic (batoid.Optic) – Optical system

  • theta_x, theta_y (float) – Field angle in radians

  • wavelength (float) – Wavelength in meters

  • nrad (int, optional) – Number of ray radii to use. (see RayVector.asPolar())

  • naz (int, optional) – Approximate number of azimuthal angles in outermost ring. (see RayVector.asPolar())

  • projection ({‘postel’, ‘zemax’, ‘gnomonic’, ‘stereographic’, ‘lambert’, ‘orthographic’}) – Projection used to convert field angle to direction cosines.

Returns:

fL (float) – Focal length in meters.

Notes

Computed from Jacobian drdth, suitably averaged over azimuth.