Rays

The most fundamental type in batoid is the RayVector. This class represents a collection of geometric rays. These rays can be generated, traced through optical systems, and then interrogated to provide insight into PSFs, aberrations, vignetting, and more.

class batoid.RayVector(x, y, z, vx, vy, vz, t=0.0, wavelength=0.0, flux=1.0, vignetted=False, failed=False, coordSys=CoordSys(array([0., 0., 0.]), array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]])))[source]

Create RayVector from 1d parameter arrays. Always makes a copy of input arrays.

Parameters:
  • x, y, z (ndarray of float, shape (n,)) – Positions of rays in meters.

  • vx, vy, vz (ndarray of float, shape (n,)) – Velocities of rays in units of the speed of light in vacuum.

  • t (ndarray of float, shape (n,)) – Reference times (divided by the speed of light in vacuum) in units of meters.

  • wavelength (ndarray of float, shape (n,)) – Vacuum wavelengths in meters.

  • flux (ndarray of float, shape (n,)) – Fluxes in arbitrary units.

  • vignetted (ndarray of bool, shape (n,)) – True where rays have been vignetted.

  • coordSys (CoordSys) – Coordinate system in which this ray is expressed. Default: the global coordinate system.

positionAtTime(t)[source]

Calculate the positions of the rays at a given time.

Parameters:

t (float) – Time (over vacuum speed of light; in meters).

Returns:

ndarray of float, shape (n, 3) – Positions in meters.

propagate(t)[source]

Propagate this RayVector to given time.

Parameters:

t (float) – Time (over vacuum speed of light; in meters).

Returns:

RayVector – Reference to self, no copy is made.

phase(r, t)[source]

Calculate plane wave phases at given position and time.

Parameters:
  • r (ndarray of float, shape (3,)) – Position in meters at which to compute phase

  • t (float) – Time (over vacuum speed of light; in meters).

Returns:

ndarray of float, shape(n,)

amplitude(r, t)[source]

Calculate (scalar) complex electric-field amplitudes at given position and time.

Parameters:
  • r (ndarray of float, shape (3,)) – Position in meters.

  • t (float) – Time (over vacuum speed of light; in meters).

Returns:

ndarray of complex, shape (n,)

sumAmplitude(r, t, ignoreVignetted=True)[source]

Calculate the sum of (scalar) complex electric-field amplitudes of all rays at given position and time.

Parameters:
  • r (ndarray of float, shape (3,)) – Position in meters.

  • t (float) – Time (over vacuum speed of light; in meters).

Returns:

complex

classmethod asGrid(optic=None, backDist=None, medium=None, stopSurface=None, coordSys=None, wavelength=None, source=None, dirCos=None, theta_x=None, theta_y=None, projection='postel', nx=None, ny=None, dx=None, dy=None, lx=None, ly=None, flux=1, nrandom=None, rng=None)[source]

Create RayVector on a parallelogram shaped region.

This function will often be used to create a grid of rays on a square grid, but is flexible enough to also create grids on an arbitrary parallelogram, or even randomly distributed across an arbitrary parallelogram-shaped region.

The algorithm starts by placing rays on the “stop” surface, and then backing them up such that they are in front of any surfaces of the optic they’re intended to trace.

The stop surface of most large telescopes is the plane perpendicular to the optic axis and flush with the rim of the primary mirror. This plane is usually also the entrance pupil since there are no earlier refractive or reflective surfaces. However, since this plane is a bit difficult to locate automatically, the default stop surface in batoid is the global x-y plane.

If a telescope has a stopSurface attribute in its yaml file, then this is usually a good choice to use in this function. Using a curved surface for the stop surface is allowed, but is usually a bad idea as this may lead to a non-uniformly illuminated pupil and is inconsistent with, say, an incoming uniform spherical wave or uniform plane wave.

Parameters:
  • optic (batoid.Optic, optional) – If present, then try to extract values for backDist, medium, stopSurface, and lx from the Optic. Note that values explicitly passed to asGrid as keyword arguments override those extracted from optic.

  • backDist (float, optional) – Map rays backwards from the stop surface to the plane that is perpendicular to the rays and backDist meters from the point (0, 0, z(0,0)) on the stop surface. This should generally be set large enough that any obscurations or phantom surfaces occuring before the stop surface are now “in front” of the ray. If this keyword is set to None and the optic keyword is set, then infer a value from optic.backDist. If both this keyword and optic are None, then use a default of 40 meters, which should be sufficiently large for foreseeable telescopes.

  • medium (batoid.Medium, optional) – Initial medium of each ray. If this keyword is set to None and the optic keyword is set, then infer a value from optic.inMedium. If both this keyword and optic are None, then use a default of vacuum.

  • stopSurface (batoid.Interface, optional) – Surface defining the system stop. If this keyword is set to None and the optic keyword is set, then infer a value from optic.stopSurface. If both this keyword and optic are None, then use a default Interface(Plane()), which is the global x-y plane.

  • wavelength (float) – Vacuum wavelength of rays in meters.

  • source (None or ndarray of float, shape (3,), optional) – Where rays originate. If None, then rays originate an infinite distance away, in which case the dirCos kwarg must also be specified to set the direction of ray propagation. If an ndarray, then the rays originate from this point in global coordinates and the dirCos kwarg is ignored.

  • dirCos (ndarray of float, shape (3,), optional) – If source is None, then this indicates the initial direction of propagation of the rays. If source is not None, then this is ignored. Also see theta_x, theta_y as an alternative to this keyword.

  • theta_x, theta_y (float, optional) – Field angle in radians. If source is None, then this indicates the initial direction of propagation of the rays. If source is not None, then this is ignored. Uses utils.fieldToDirCos to convert to direction cosines. Also see dirCos as an alternative to this keyword.

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

  • nx, ny (int, optional) – Number of rays on each side of grid.

  • dx, dy (float or (2,) array of float, optional) – Separation in meters between adjacent rays in grid. If scalars, then the separations are exactly along the x and y directions. If arrays, then these are interpretted as the primitive vectors for the first and second dimensions of the grid. If only dx is explicitly specified, then dy will be inferred as a 90-degree rotation from dx with the same length as dx.

  • lx, ly (float or (2,) array of float, optional) – Length of each side of ray grid. If scalars, then these are measured along the x and y directions. If arrays, then these also indicate the primitive vectors orientation of the grid. If only lx is specified, then ly will be inferred as a 90-degree rotation from lx with the same length as lx. If lx is None, then first infer a value from nx and dx, and if that doesn’t work, infer a value from optic.pupilSize.

  • flux (float, optional) – Flux to assign each ray. Default is 1.0.

  • nrandom (None or int, optional) – If not None, then uniformly sample this many rays from parallelogram region instead of sampling on a regular grid.

  • rng (None or int or numpy.random.Generator, optional) – Random number generator or seed to use for random sampling.

classmethod asPolar(optic=None, backDist=None, medium=None, stopSurface=None, coordSys=None, wavelength=None, outer=None, inner=None, source=None, dirCos=None, theta_x=None, theta_y=None, projection='postel', nrad=None, naz=None, flux=1, nrandom=None, rng=None)[source]

Create RayVector on an annular region using a hexapolar grid.

This function can be used to regularly sample the entrance pupil of a telescope using polar symmetry (really, hexagonal symmetry). Rings of different radii are used, with the number of samples on each ring restricted to a multiple of 6 (with the exception of a potential central “ring” of radius 0, which is only ever sampled once). This may be more efficient than using a square grid since more of the rays generated may avoid vignetting.

This function is also used to generate rays uniformly randomly sampled from a given annular region.

The algorithm used here starts by placing rays on the “stop” surface, and then backing them up such that they are in front of any surfaces of the optic they’re intended to trace.

The stop surface of most large telescopes is the plane perpendicular to the optic axis and flush with the rim of the primary mirror. This plane is usually also the entrance pupil since there are no earlier refractive or reflective surfaces. However, since this plane is a bit difficult to locate automatically, the default stop surface in batoid is the global x-y plane.

If a telescope has a stopSurface attribute in its yaml file, then this is usually a good choice to use in this function. Using a curved surface for the stop surface is allowed, but is usually a bad idea as this may lead to a non-uniformly illuminated pupil and is inconsistent with, say, an incoming uniform spherical wave or uniform plane wave.

Parameters:
  • optic (batoid.Optic, optional) – If present, then try to extract values for backDist, medium, stopSurface, and outer from the Optic. Note that values explicitly passed to asPolar as keyword arguments override those extracted from optic.

  • backDist (float, optional) – Map rays backwards from the stop surface to the plane that is perpendicular to the ray and backDist meters from the point (0, 0, z(0,0)) on the stop surface. This should generally be set large enough that any obscurations or phantom surfaces occuring before the stop surface are now “in front” of the ray. If this keyword is set to None and the optic keyword is set, then infer a value from optic.backDist. If both this keyword and optic are None, then use a default of 40 meters, which should be sufficiently large for foreseeable telescopes.

  • medium (batoid.Medium, optional) – Initial medium of each ray. If this keyword is set to None and the optic keyword is set, then infer a value from optic.inMedium. If both this keyword and optic are None, then use a default of vacuum.

  • stopSurface (batoid.Interface, optional) – Surface defining the system stop. If this keyword is set to None and the optic keyword is set, then infer a value from optic.stopSurface. If both this keyword and optic are None, then use a default Interface(Plane()), which is the global x-y plane.

  • wavelength (float) – Vacuum wavelength of rays in meters.

  • outer (float) – Outer radius of annulus in meters.

  • inner (float, optional) – Inner radius of annulus in meters. Default is 0.0.

  • source (None or ndarray of float, shape (3,), optional) – Where rays originate. If None, then rays originate an infinite distance away, in which case the dirCos kwarg must also be specified to set the direction of ray propagation. If an ndarray, then the rays originate from this point in global coordinates and the dirCos kwarg is ignored.

  • dirCos (ndarray of float, shape (3,), optional) – If source is None, then this indicates the initial direction of propagation of the rays. If source is not None, then this is ignored. Also see theta_x, theta_y as an alternative to this keyword.

  • theta_x, theta_y (float, optional) – Field angle in radians. If source is None, then this indicates the initial direction of propagation of the rays. If source is not None, then this is ignored. Uses utils.fieldToDirCos to convert to direction cosines. Also see dirCos as an alternative to this keyword.

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

  • nrad (int) – Number of radii on which create rays.

  • naz (int) – Approximate number of azimuthal angles uniformly spaced along the outermost ring. Each ring is constrained to have a multiple of 6 azimuths, so the realized value may be slightly different than the input value here. Inner rings will have fewer azimuths in proportion to their radius, but will still be constrained to a multiple of 6. (If the innermost ring has radius 0, then exactly 1 ray, with azimuth undefined, will be used on that “ring”.)

  • flux (float, optional) – Flux to assign each ray. Default is 1.0.

  • nrandom (int, optional) – If not None, then uniformly sample this many rays from annular region instead of sampling on a hexapolar grid.

  • rng (None or int or numpy.random.Generator, optional) – Random number generator or seed to use for random sampling.

classmethod asSpokes(optic=None, backDist=None, medium=None, stopSurface=None, coordSys=None, wavelength=None, outer=None, inner=0.0, source=None, dirCos=None, theta_x=None, theta_y=None, projection='postel', spokes=None, rings=None, spacing='uniform', flux=1)[source]

Create RayVector on an annular region using a spokes pattern.

The function generates rays on a rings-and-spokes pattern, with a fixed number of radii for each azimuth and a fixed number of azimuths for each radius. Its main use is for decomposing functions in pupil space into Zernike components using Gaussian Quadrature integration on annuli. For more general purpose annular sampling, RayVector.asPolar() is often a better choice since it samples the pupil more uniformly.

The algorithm used here starts by placing rays on the “stop” surface, and then backing them up such that they are in front of any surfaces of the optic they’re intended to trace.

The stop surface of most large telescopes is the plane perpendicular to the optic axis and flush with the rim of the primary mirror. This plane is usually also the entrance pupil since there are no earlier refractive or reflective surfaces. However, since this plane is a bit difficult to locate automatically, the default stop surface in batoid is the global x-y plane.

If a telescope has a stopSurface attribute in its yaml file, then this is usually a good choice to use in this function. Using a curved surface for the stop surface is allowed, but is usually a bad idea as this may lead to a non-uniformly illuminated pupil and is inconsistent with, say, an incoming uniform spherical wave or uniform plane wave.

Parameters:
  • optic (batoid.Optic, optional) – If present, then try to extract values for backDist, medium, stopSurface, and outer from the Optic. Note that values explicitly passed to asSpokes as keyword arguments override those extracted from optic.

  • backDist (float, optional) – Map rays backwards from the stop surface to the plane that is perpendicular to the ray and backDist meters from the point (0, 0, z(0,0)) on the stop surface. This should generally be set large enough that any obscurations or phantom surfaces occuring before the stop surface are now “in front” of the ray. If this keyword is set to None and the optic keyword is set, then infer a value from optic.backDist. If both this keyword and optic are None, then use a default of 40 meters, which should be sufficiently large for foreseeable telescopes.

  • medium (batoid.Medium, optional) – Initial medium of each ray. If this keyword is set to None and the optic keyword is set, then infer a value from optic.inMedium. If both this keyword and optic are None, then use a default of vacuum.

  • stopSurface (batoid.Interface, optional) – Surface defining the system stop. If this keyword is set to None and the optic keyword is set, then infer a value from optic.stopSurface. If both this keyword and optic are None, then use a default Interface(Plane()), which is the global x-y plane.

  • wavelength (float) – Vacuum wavelength of rays in meters.

  • outer (float) – Outer radius of annulus in meters.

  • inner (float, optional) – Inner radius of annulus in meters. Default is 0.0.

  • source (None or ndarray of float, shape (3,), optional) – Where rays originate. If None, then rays originate an infinite distance away, in which case the dirCos kwarg must also be specified to set the direction of ray propagation. If an ndarray, then the rays originate from this point in global coordinates and the dirCos kwarg is ignored.

  • dirCos (ndarray of float, shape (3,), optional) – If source is None, then this indicates the initial direction of propagation of the rays. If source is not None, then this is ignored. Also see theta_x, theta_y as an alternative to this keyword.

  • theta_x, theta_y (float, optional) – Field angle in radians. If source is None, then this indicates the initial direction of propagation of the rays. If source is not None, then this is ignored. Uses utils.fieldToDirCos to convert to direction cosines. Also see dirCos as an alternative to this keyword.

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

  • spokes (int or ndarray of float) – If int, then number of spokes to use. If ndarray, then the values of the spokes azimuthal angles in radians.

  • rings (int or ndarray of float) – If int, then number of rings to use. If array, then the values of the ring radii to use in meters.

  • spacing ({‘uniform’, ‘GQ’}) – If uniform, assign ring radii uniformly between inner and outer. If GQ, then assign ring radii as the Gaussian Quadrature points for integration on an annulus. In this case, the ray fluxes will be set to the Gaussian Quadrature weights (and the flux kwarg will be ignored).

  • flux (float, optional) – Flux to assign each ray. Default is 1.0.

classmethod fromStop(x, y, optic=None, backDist=None, medium=None, stopSurface=None, coordSys=None, wavelength=None, source=None, dirCos=None, theta_x=None, theta_y=None, projection='postel', flux=1)[source]

Create rays that intersects the “stop” surface at given points.

The algorithm used here starts by placing the rays on the “stop” surface, and then backing them up such that they are in front of any surfaces of the optic they’re intended to trace.

The stop surface of most large telescopes is the plane perpendicular to the optic axis and flush with the rim of the primary mirror. This plane is usually also the entrance pupil since there are no earlier refractive or reflective surfaces. However, since this plane is a bit difficult to locate automatically, the default stop surface in batoid is the global x-y plane.

If a telescope has a stopSurface attribute in its yaml file, then this is usually a good choice to use in this function. Using a curved surface for the stop surface is allowed, but is usually a bad idea as this may lead to a non-uniformly illuminated pupil and is inconsistent with, say, an incoming uniform spherical wave or uniform plane wave.

Parameters:
  • x, y (ndarray) – X/Y coordinates on the stop surface where the rays would intersect if not refracted or reflected first.

  • optic (batoid.Optic, optional) – If present, then try to extract values for backDist, medium, and stopSurface from the Optic. Note that values explicitly passed here as keyword arguments override those extracted from optic.

  • backDist (float, optional) – Map rays backwards from the stop surface to the plane that is perpendicular to the rays and backDist meters from the point (0, 0, z(0,0)) on the stop surface. This should generally be set large enough that any obscurations or phantom surfaces occuring before the stop surface are now “in front” of the ray. If this keyword is set to None and the optic keyword is set, then infer a value from optic.backDist. If both this keyword and optic are None, then use a default of 40 meters, which should be sufficiently large for foreseeable telescopes.

  • medium (batoid.Medium, optional) – Initial medium of rays. If this keyword is set to None and the optic keyword is set, then infer a value from optic.inMedium. If both this keyword and optic are None, then use a default of vacuum.

  • stopSurface (batoid.Interface, optional) – Surface defining the system stop. If this keyword is set to None and the optic keyword is set, then infer a value from optic.stopSurface. If both this keyword and optic are None, then use a default Interface(Plane()), which is the global x-y plane.

  • wavelength (float) – Vacuum wavelength of rays in meters.

  • source (None or ndarray of float, shape (3,), optional) – Where the rays originate. If None, then the rays originate an infinite distance away, in which case the dirCos kwarg must also be specified to set the direction of ray propagation. If an ndarray, then the rays originates from this point in global coordinates and the dirCos kwarg is ignored.

  • dirCos (ndarray of float, shape (3,), optional) – If source is None, then indicates the direction of ray propagation. If source is not None, then this is ignored.

  • theta_x, theta_y (float, optional) – Field angle in radians. If source is None, then this indicates the initial direction of propagation of the rays. If source is not None, then this is ignored. Uses utils.fieldToDirCos to convert to direction cosines. Also see dirCos as an alternative to this keyword.

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

  • flux (float, optional) – Flux of rays. Default is 1.0.

classmethod fromFieldAngles(theta_x, theta_y, projection='postel', optic=None, backDist=None, medium=None, stopSurface=None, wavelength=None, x=0, y=0, flux=1)[source]

Create RayVector with one stop surface point but many field angles.

This method is similar to fromStop but broadcasts over theta_x and theta_y instead of over x and y. There is less currently less effort paid to synchronizing the t values of the created rays, as they don’t correspond to points on a physical incoming wavefront in this case. The primary intended use case is to map chief rays (x = y = 0) from incoming field angle to focal plane position.

Parameters:
  • theta_x, theta_y (ndarray) – Field angles in radians.

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

  • optic (batoid.Optic, optional) – If present, then try to extract values for backDist, medium, and stopSurface from the Optic. Note that values explicitly passed here as keyword arguments override those extracted from optic.

  • backDist (float, optional) – Map rays backwards from the stop surface this far. This should generally be set large enough that any obscurations or phantom surfaces occuring before the stop surface are now “in front” of the rays. If this keyword is set to None and the optic keyword is set, then infer a value from optic.backDist. If both this keyword and optic are None, then use a default of 40 meters, which should be sufficiently large for foreseeable telescopes.

  • medium (batoid.Medium, optional) – Initial medium of rays. If this keyword is set to None and the optic keyword is set, then infer a value from optic.inMedium. If both this keyword and optic are None, then use a default of vacuum.

  • stopSurface (batoid.Interface, optional) – Surface defining the system stop. If this keyword is set to None and the optic keyword is set, then infer a value from optic.stopSurface. If both this keyword and optic are None, then use a default Interface(Plane()), which is the global x-y plane.

  • wavelength (float) – Vacuum wavelength of rays in meters.

  • x, y (float) – X/Y coordinates on the stop surface where the rays would intersect if not refracted or reflected first.

  • flux (float, optional) – Flux of rays. Default is 1.0.

property r

Positions of rays in meters.

Type:

ndarray of float, shape (n, 3)

property x

The x components of ray positions in meters.

property y

The y components of ray positions in meters.

property z

The z components of ray positions in meters.

property v

Velocities of rays in units of the speed of light in vacuum. Note that these may have magnitudes < 1 if the rays are inside a refractive medium.

Type:

ndarray of float, shape (n, 3)

property vx

The x components of ray velocities units of the vacuum speed of light.

property vy

The y components of ray velocities units of the vacuum speed of light.

property vz

The z components of ray velocities units of the vacuum speed of light.

property t

Reference times (divided by the speed of light in vacuum) in units of meters, also known as the optical path lengths.

property wavelength

Vacuum wavelengths in meters.

property flux

Fluxes in arbitrary units.

property vignetted

True for rays that have been vignetted.

property failed

True for rays that have failed. This may occur, for example, if batoid failed to find the intersection of a ray wiht a surface.

property k

Wavevectors of plane waves in units of radians per meter. The magnitude of each wavevector is equal to \(2 \pi n / \lambda\), where \(n\) is the refractive index and \(\lambda\) is the wavelength.

Type:

ndarray of float, shape (n, 3)

property kx

The x component of each ray wavevector in radians per meter.

property ky

The y component of each ray wavevector in radians per meter.

property kz

The z component of each ray wavevector in radians per meter.

property omega

The temporal angular frequency of each plane wave divided by the vacuum speed of light in units of radians per meter. Equals \(2 \pi / \lambda\).

toCoordSys(coordSys)[source]

Transform this RayVector into a new coordinate system.

Parameters:

coordSys (batoid.CoordSys) – Destination coordinate system.

Returns:

RayVector – Reference to self, no copy is made.