psfsim.romantrace
Attributes
Classes
Class defining a ray bundle, constructed from a field position. |
Functions
|
Lanczos function. |
|
Apply Lanczos reweighting to boundary cells using high-resolution data. |
|
Function to compute the index of refraction of Infrasil 301. |
|
Makes a transformation matrix from a set of translations and rotations. |
|
Carries out trace through RST optics. |
|
Carries out trace through RST optics. |
|
Demo and test functions for romantrace. |
Module Contents
- _lanczos_weight(dx, dy, a=3)[source]
Lanczos function.
- Parameters:
dx (float or np.ndarray) – The input x value(s).
dy (float or np.ndarray,) – The input y value(s). Must be the same shape as dx.
a (int, optional) – The Lanczos parameter (default is 3).
- Returns:
The Lanczos function evaluated on the dx by dy grid.
- Return type:
float or np.ndarray
- _apply_lanczos_reweighting(RB_open_lowres, RB_open_hires, bdycells, ovsamp, a_lanczos)[source]
Apply Lanczos reweighting to boundary cells using high-resolution data.
- Parameters:
RB_open_lowres (ndarray) – Low-resolution aperture array
RB_open_hires (ndarray) – High-resolution aperture array (flattened to num_bdy, ovsamp, ovsamp)
bdycells (tuple) – (x_indices, y_indices) of boundary cells
ovsamp (int) – Oversampling factor
a_lanczos (int) – Lanczos kernel size parameter
- Returns:
new_values – Reweighted aperture values at boundary cells
- Return type:
ndarray
- n_Infrasil301(wl, T=180.0)[source]
Function to compute the index of refraction of Infrasil 301.
- Parameters:
wl (float) – Vacuum wavelength in millimeters.
T (float, optional) – Temperature in Kelvin.
- Returns:
The (real) index of refraction.
- Return type:
float
Notes
Original reference is: Sellmeier coefficients from Leviton, Frey, & Madison (2008) arXiv:0805.0096
Valid in the Roman wavelength range.
- build_transform_matrix(xde=0.0, yde=0.0, zde=0.0, ade=0.0, bde=0.0, cde=0.0, unit='degree')[source]
Makes a transformation matrix from a set of translations and rotations.
The matrix is 4x4 and the convention is:
u(global coords) = R u(obj coords) where u is a vector of [1,x,y,z]
- Parameters:
xde (float, optional) – Translations on each axis. The convention is that the origin in object coordinates is at (xde, yde, zde) in global coordinates.
yde (float, optional) – Translations on each axis. The convention is that the origin in object coordinates is at (xde, yde, zde) in global coordinates.
zde (float, optional) – Translations on each axis. The convention is that the origin in object coordinates is at (xde, yde, zde) in global coordinates.
ade (float, optional) –
Euler angles, in the rotX - rotY - rotZ convention (so you can envision rotating the global to the object system via a right-handed rotation around X; then a right-handed rotation around Y; and finally a left-handed rotation around Z). The convention is that:
The object z-axis direction is in the global direction (-sin bde, sin ade * cos bde, cos ade * cos bde)
The global x-axis direction is in the object direction (cos bde * cos cde, cos bde * sin cde, -sin bde).
bde (float, optional) –
Euler angles, in the rotX - rotY - rotZ convention (so you can envision rotating the global to the object system via a right-handed rotation around X; then a right-handed rotation around Y; and finally a left-handed rotation around Z). The convention is that:
The object z-axis direction is in the global direction (-sin bde, sin ade * cos bde, cos ade * cos bde)
The global x-axis direction is in the object direction (cos bde * cos cde, cos bde * sin cde, -sin bde).
cde (float, optional) –
Euler angles, in the rotX - rotY - rotZ convention (so you can envision rotating the global to the object system via a right-handed rotation around X; then a right-handed rotation around Y; and finally a left-handed rotation around Z). The convention is that:
The object z-axis direction is in the global direction (-sin bde, sin ade * cos bde, cos ade * cos bde)
The global x-axis direction is in the object direction (cos bde * cos cde, cos bde * sin cde, -sin bde).
unit (str, optional) – Angular unit: “degree” (default) or “radian”.
- Returns:
The 4x4 transformation matrix.
- Return type:
np.ndarray
- class RayBundle(xan, yan, N, hasE=False, width=2500.0, startpos=3500.0, wl=0.00129, wlref=0.00129, jacobian=None, hires=None, ovsamp=6, idealgeom=True)[source]
Class defining a ray bundle, constructed from a field position.
- Parameters:
xan (float) – Field angles (in degrees).
yan (float) – Field angles (in degrees).
N (int) – The bundle size on each axis of the entrance pupil (so
N**2total points).hasE (bool, optional) – Also propagate the electric field?
width (float, optional) – The size of the grid for the entrance pupil (in mm; default of 2500 mm is good for Roman).
startpos (float, optional) – Starting Z-position in mm (default is before the SM obstruction in Roman).
wl (float, optional) – The vacuum wavelength in mm.
wlref (float, optional) – The vacuum wavelength in mm used as a reference (this is so that when we make chromatic PSFs, the “center” of the PSF always ends up in the same place, and e.g. DCR effects show up in wavelength-dependent Zernike Z2 & Z3 modes).
jacobian (np.ndarray, optional) – If used, this will give a 2x2 distortion matrix, used so that the output exit pupil is on a square grid. Default is a square grid on the entrance pupil.
hires (list of np.ndarray of int, optional) – If given, then instead of the full grid, trace only the cells given;
hires[0]is a 1D array of y-values andhires[1]is a 1D array of x-values.ovsamp (int, optional) – Oversamples cells in the entrance pupil by this factor; only used if hires is given.
idealgeom (bool, optional) – Forces the design model rather than with the best-fit offsets.
- N1, N2
Trimmed size (needed for hires mode; otherwise equal to N).
- Type:
int
- x[source]
3D array of position of rays, shape (N1, N2, 4); unit = mm; the last axis is a position so has
x[:,:,0] == 1.- Type:
np.ndarray of float
- p[source]
3D array of direction of rays, shape (N1, N2, 4); the last axis is a direction so has
p[:,:,0] == 0.- Type:
np.ndarray of float
- xan, yan
Field angles (in degrees).
- Type:
float
- wl, wlref
Wavelength and reference wavelength (the latter for geometric trace). Units are mm.
- Type:
float
- E
The electric field (optional, 4D: 2D for array, 1D for input pol, 1D for output pol). Shape is (N1, N2, 2, 2). None if not used.
- Type:
np.ndarray of complex or None.
- classmethod MV(a, b)[source]
Matrix-vector multiplication.
- Parameters:
a (np.ndarray) – Matrix, shape (4, 4).
b (np.ndarray) – Array of vectors, shape (…, 4).
- Returns:
Product a times b, acting on last axis, same shape as b.
- Return type:
np.ndarray
See also
MiVThis is the inverse function.
- classmethod MiV(a, b)[source]
Inverse matrix-vector multiplication.
- Parameters:
a (np.ndarray) – Matrix, shape (4, 4).
b (np.ndarray) – Array of vectors, shape (…, 4).
- Returns:
Product inv(a) times b, acting on last axis, same shape as b.
- Return type:
np.ndarray
See also
MVThe forward function (
MiVis the inverse).
- intersect_surface(Trf, Rinv=0.0, K=0.0, update=True)[source]
Gets intersection of a ray bundle and a conic section surface.
Updates the path with intersection information (unless update is set to False, in which case just the intersection geometry is returned but the ray positions don’t update to the intersection surface).
- Parameters:
Trf (np.ndarray of float) – The 4x4 transformation matrix for locating the surface.
Rinv (float, optional) – The inverse radius of curvature of the surface; 0 = flat (default), positive = curves toward object +Z, negative = curves toward object -Z.
K (float, optional) – Conic constant (-K = eccentricity^2), 0 = spherical (default), -1 = parabola.
update (bool, optional) – If set to False, only check where the rays hit the surface.
- Returns:
pos_object (np.ndarray of float) – Object coordinates of where the rays hit the surface; shape (N, N, 2) (so only object x and y coordinates are returned).
dir_global (np.ndarray of float) – Global unit normal vector of the surface; shape (N, N, 4). Since this is a direction, ``dir_global[:, :, 0] == 0`, and then the 1, 2, and 3 components correspond to the x, y, and z directions.
L (np.ndarray of float) – 2D array of the path lengths traversed to reach the surface.
- mask(Trf, Rinv, K, R, masklist)[source]
Masks incoming rays at a given surface.
- Parameters:
Trf (np.ndarray of float) – The 4x4 transformation matrix for the surface.
Rinv (float) – The inverse radius of curvature of the surface.
K (float) – The conic constant for the surface.
R (float or None) – If not None, mask rays whose intersection points are outside radius R.
masklist (list of dict) –
A list of obstructions on the surface. The attributes are:
CIR,REX,REY(CODE V codes for shapes)ADX,ADY,ARO(CODE V de-centers)
- Return type:
None
Notes
All dimensions are in millimeters.
- intersect_surface_and_reflect(Trf, Rinv=0.0, K=0.0, rCoefs=None, activeZone=None, gd=None)[source]
Propagates rays to a surface and performs a reflection.
Paramters
- Trfnp.ndarray of float
The 4x4 transformation matrix for the surface.
- Rinvfloat, optional
The inverse radius of curvature of the surface. Default = plane.
- Kfloat, optional
The conic constant for the surface. Default = sphere.
- rCoefsfunction, optional
If given, this is a function that takes in an array of incidence angles (1st argument) and a vacuum wavelength (2nd argument) and returns S- and P-polarized reflection coefficient arrays (complex: same shape as angle of incidence). Otherwise assumes a perfect reflecting condition.
- activeZonedict, optional
A dictionary of CODEV codes for the active region. Valid keys are CIR`,
REX,REY,ADX,ADY, andARO.- gddict, optional
A dictionary of parameters for gradient computation. Should have keys:
grad: where to write the gradient map; andsurface: which surface (e.g.,M1).
- rtype:
None
See also
intersect_surfaceThis function differs in that it also performs the reflection.
- intersect_surface_and_refract(Trf, Rinv=0.0, K=0.0, n_new=1.0, tCoefs=None, activeZone=None, gd=None)[source]
Propagates rays to a surface and performs a refraction.
Paramters
- Trfnp.ndarray of float
The 4x4 transformation matrix for the surface.
- Rinvfloat, optional
The inverse radius of curvature of the surface. Default = plane.
- Kfloat, optional
The conic constant for the surface. Default = sphere.
- n_newfloat, optional
The index of refraction of the new medium. Default = vacuum.
- tCoefsfunction, optional
If given, this is a function that takes in an array of incidence angles (1st argument) and a vacuum wavelength (2nd argument) and returns S- and P-polarized reflection coefficient arrays (complex: same shape as angle of incidence). Otherwise assumes a perfect reflecting condition.
- activeZonedict, optional
A dictionary of CODEV codes for the active region. Valid keys are CIR`,
REX,REY,ADX,ADY, andARO.- gddict, optional
A dictionary of parameters for gradient computation. Should have keys:
grad: where to write the gradient map; andsurface: which surface (e.g.,M1).
- rtype:
None
See also
intersect_surfaceThis function differs in that it also performs the refraction.
- _RomanRayBundle(xan, yan, N, usefilter, wl=None, hasE=False, width=2500.0, jacobian=None, hires=None, ovsamp=6, idealgeom=False, idealmirror=False, outsca=None, errs=None)[source]
Carries out trace through RST optics.
- Parameters:
xan (float) – Angles in degrees in WFI local field angles.
yan (float) – Angles in degrees in WFI local field angles.
N (int) – Pupil sampling (NxN grid).
usefilter (char) – The RST filter element. One of ‘R’, ‘Z’, ‘Y’, ‘J’, ‘H’, ‘F’, ‘K’, ‘W’.
wl (float, optional) – The wavelength in mm. Default is the central wavelength of that filter.
hasE (bool, optional) – If True, propagate the electric field.
width (float, optional) – The size of the grid for sampling of the entrance pupil, in mm. (Default is good for RST.)
jacobian (np.ndarray, optional) – If used, this will give a 2x2 distortion matrix, used so that the output exit pupil is on a square grid. Default is a square grid on the entrance pupil.
hires (list of np.ndarray of int, optional) – If given, then instead of the full grid, trace only the cells given;
hires[0]is a 1D array of y-values andhires[1]is a 1D array of x-values.ovsamp (int, optional) – Oversamples cells in the entrance pupil by this factor; only used if hires is given.
idealgeom (bool, optional) – Forces the design model rather than with the best-fit offsets.
idealmirror (bool, optional) – Forces the mirror to be an ideal conducting surface instead of the model.
outsca (int, optional) – Which output SCA to project the PSF on to? (1..18) The default is to choose the SCA whose center is closest to where the ray lands.
errs (dict, optional) –
The surface error model. Dictionary with keys:
basis:
RomanBasisSetobject.grad : bool, optional: whether to compute gradient of s with respect to the surface error parameters.
Warning: this can be big!
arr: np.ndarray, optional: 1D array of amplitudes of each basis mode.
- Returns:
The ray bundle object. The following additional information is provided as numpy arrays:
RB.x_out : shape (2,), location of output ray on FPA [in mm] RB.xyi : shape (N,N,2), coordinates of the initial rays [in mm] RB.u : shape (N,N,2), directions (orthographic projection) RB.s : shape (N,N), optical path length of ray to position s
If hasE is True, then also provides:
RB.E : shape (N,N,2,4), complex, electric field for the 2 initial polarizations and 3 components
(last axis 0th component should be 0)
- If
errs["grad"]is True, then also provides: RB.grad : shape (N,N,basis_set.N), derivative of s with respect to each surface error parameter.
- If
- Return type:
- RomanRayBundle(xan, yan, N, usefilter, wl=None, hasE=False, width=2500.0, jacobian=None, ovsamp=6, a_lanczos=3, idealgeom=False, idealmirror=False, outsca=None, errs=None)[source]
Carries out trace through RST optics.
- Parameters:
xan (float) – Angles in degrees in WFI local field angles.
yan (float) – Angles in degrees in WFI local field angles.
N (int) – Pupil sampling (NxN grid).
usefilter (char) – The RST filter element. One of ‘R’, ‘Z’, ‘Y’, ‘J’, ‘H’, ‘F’, ‘K’, ‘W’.
wl (float, optional) – The wavelength in mm. Default is the central wavelength of that filter.
hasE (bool, optional) – If True, propagate the electric field.
width (float, optional) – The size of the grid for sampling of the entrance pupil, in mm. (Default is good for RST.)
jacobian (np.ndarray, optional) – If used, this will give a 2x2 distortion matrix, used so that the output exit pupil is on a square grid. Default is a square grid on the entrance pupil.
ovsamp (int, optional) – Oversamples cells in the entrance pupil by this factor; only used if hires is given.
a_lanczos (int, optional) – The “a” parameter for Lanczos interpolation; this controls the size of the kernel. Default is 3.
idealgeom (bool, optional) – Forces the design model rather than with the best-fit offsets.
idealmirror (bool, optional) – Forces the mirror to be an ideal conducting surface instead of the model.
outsca (int, optional) – Which output SCA to project the PSF on to? (1..18) The default is to choose the SCA whose center is closest to where the ray lands.
errs (dict, optional) –
The surface error model. Dictionary with keys:
basis:
RomanBasisSetobject.grad : bool, optional: whether to compute gradient of s with respect to the surface error parameters.
Warning: this can be big!
arr: np.ndarray, optional: 1D array of amplitudes of each basis mode.
- Returns:
The ray bundle object. The following additional information is provided as numpy arrays:
RB.x_out : shape (2,), location of output ray on FPA [in mm] RB.xyi : shape (N,N,2), coordinates of the initial rays [in mm] RB.u : shape (N,N,2), directions (orthographic projection) RB.s : shape (N,N), optical path length of ray to position s RS.open : shape (N,N), fraction of that cell that is open (between 0 and 1 inclusive).
If hasE is True, then also provides:
RB.E : shape (N,N,2,4), complex, electric field for the 2 initial polarizations and 3 components
(last axis 0th component should be 0)
- If
errs["grad"]is True, then also provides: RB.grad : shape (N,N,basis_set.N), derivative of s with respect to each surface error parameter.
- If
- Return type:
See also
_RomanRayBundleThe routine “wrapped inside” (this is called once to figure out which cells need to be simulated at higher resolution).
Notes
This is based on the specifications in: “Opto-Mechanical Definitions”, RST-SYS-SPEC-0055, Revision E released by the Configuration Management Office June 11, 2021 (not export controlled). Most information is in CODE V format, but some was converted to be usable in this Python script.