Documentation for pysiaf¶
SIAF files contain detailed geometric focal plane descriptions and relationships for the science instruments. They are maintained in the JWST/HST PRD (Project Reference Database). pysiaf is a python package to access, interpret, maintain, and generate SIAF, in particular for JWST. Tools for applying the frame transformations, plotting, comparison, and validation are provided.
pysiaf¶
Handling of Science Instrument Aperture Files (SIAF) for space telescopes
Functionalities¶
Captures current PRD content, i.e. pysiaf includes a copy of the SIAF XML files. These are maintained to be synchronized with the PRD.
Transformations between the SIAF frames (Detector, Science, Ideal, Telelescope/V) are pre-loaded and easily accessible.
Tools for plotting, validation, and comparison of SIAF apertures and files.
Support for implementing transformations between celestial (RA, Dec) and telescope/V (V2, V3) coordinate systems is provided.
Input/output: reading SIAF XML, writing XML/Excel/csv etc.
Captures SI source data and code to generate the SIAF apertures
Standard python package with installation script, unit tests, documentation.
Supports working with HST SIAF (read-only).
Where to Find pysiaf¶
pysiaf is hosted and developed at https://github.com/spacetelescope/pysiaf
Installation¶
This package is being developed in a python 3.6 environment.
How to install¶
pysiaf is available on PyPI and is included in astroconda.
pip install pysiaf
Clone the repository: git clone https://github.com/spacetelescope/pysiaf Install pysiaf: cd pysiaf python setup.py install or pip install .
Known installation issue¶
If you get an error upon import pysiaf
that traces back to import lxml.etree as ET
and states
ImportError [...] Library not loaded: libxml2.2.dylib Reason: Incompatible library version: etree.[...] requires version 12.0.0 or later, but libxml2.2.dylib provides version 10.0.0
,
this can probably be fixed by downgrading the version of lxml, e.g.
pip uninstall lxml
pip install lxml==3.6.4
User Documentation¶
Example usage:
Check which PRD version is in use:
print(pysiaf.JWST_PRD_VERSION)
Frame transformations (det
, sci
, idl
, tel
are supported frames):
import pysiaf
instrument = 'NIRISS'
# read SIAFXML
siaf = pysiaf.Siaf(instrument)
# select single aperture by name
nis_cen = siaf['NIS_CEN']
# access SIAF parameters
print('{} V2Ref = {}'.format(nis_cen.AperName, nis_cen.V2Ref))
print('{} V3Ref = {}'.format(nis_cen.AperName, nis_cen.V3Ref))
for attribute in ['InstrName', 'AperShape']:
print('{} {} = {}'.format(nis_cen.AperName, attribute, getattr(nis_cen, attribute)))
# coordinates in Science frame
sci_x = np.array([0, 2047, 2047, 0])
sci_y = np.array([0, 0, 2047, 2047])
# transform from Science frame to Ideal frame
idl_x, idl_y = nis_cen.sci_to_idl(sci_x, sci_y)
Using sky transforms¶
Transformations to/from sky
coordinates (RA, Dec) are also supported, but you have to first define and set an
attitude matrix that represents the observatory orientation with respect to the celestial sphere. This can be done with
pysiaf.utils.rotations.attitude_matrix
:
# find attitude with some coordinates (v2,v3) pointed at (ra, dec) with a given pa
attmat = pysiaf.utils.rotations.attitude_matrix(v2, v3, ra, dec, pa)
# set that attitude for the transforms
nis_cen.set_attitude_matrix(attmat)
# transform from Science frame to Sky frame
sky_ra, sky_dec = nis_cen.sci_to_sky(sci_x, sci_y)
Sky coordinates are given in units of degrees RA and Dec.
Reporting Issues / Contributing¶
Do you have feedback and feature requests? Is there something missing you would like to see? Please open a new issue or new pull request at https://github.com/spacetelescope/pysiaf for bugs, feedback, or new features you would like to see. If there is an issue you would like to work on, please leave a comment and we will be happy to assist. New contributions and contributors are very welcome! This package follows the STScI Code of Conduct strives to provide a welcoming community to all of our users and contributors.
Coding and other guidelines¶
We strive to adhere to the STScI Style Guides.
How to make a code contribution¶
The following describes the typical work flow for contributing to the pysiaf project (adapted from https://github.com/spacetelescope/jwql):
Do not commit any sensitive information (e.g. STScI-internal path structures, machine names, user names, passwords, etc.) to this public repository. Git history cannot be erased.
Create a fork off of the
spacetelescope
pysiaf
repository on your personal github space.Make a local clone of your fork.
Ensure your personal fork is pointing
upstream
to https://github.com/spacetelescope/pysiafOpen an issue on
spacetelescope
pysiaf
that describes the need for and nature of the changes you plan to make. This is not necessary for minor changes and fixes.Create a branch on that personal fork.
Make your software changes.
Push that branch to your personal GitHub repository, i.e. to
origin
.On the
spacetelescope
pysiaf
repository, create a pull request that merges the branch intospacetelescope:master
.Assign a reviewer from the team for the pull request, if you are allowed to do so.
Iterate with the reviewer over any needed changes until the reviewer accepts and merges your branch.
Delete your local copy of your branch.
Disclaimer¶
All parameter values in pysiaf are subject to change. JWST values are preliminary until the JWST observatory commissioning has concluded.
Distortion and other transformations in pysiaf are of sufficient accuracy for operations, but do not necessarily have science-grade quality. For instance, generally only one filter solution is carried per aperture. For science-grade transformations, please consult the science pipelines and their reference files (see https://jwst-docs.stsci.edu/display/JDAT/JWST+Data+Reduction+Pipeline).
For science observation planning, the focal plane geometry implemented in the latest APT (http://www.stsci.edu/hst/proposing/apt) takes precedence.
The STScI Telescopes Branch provides full support of pysiaf for S&OC operational systems only.
Reference API¶
aperture¶
Module to handle SIAF apertures.
The aperture module defined classes and functions to support working with SIAF apertures. The main class is Aperture, JwstAperture and HstAperture inherit from it. The class methods support transformations between any of the 4 SIAF frames (detector, science/DMS, ideal, V-frame/telescope). Aperture attributes correspond to entries in the SIAF xml files in the Project Reference Database (PRD), and are checked for format compliance. Methods for aperture validation and plotting are provided. The nomenclature is adapted from the JWST SIAF document Cox & Lallo: JWST-STScI-001550.
Notes¶
Numerous contributions and code snippets by Colin Cox were incorporated. Some methods were adapted from the jwxml package (https://github.com/mperrin/jwxml). Some of the polynomial transformation code was adapted from (https://github.com/spacetelescope/ramp_simulator/).
- class pysiaf.aperture.Aperture[source]¶
A class for aperture definitions of astronomical telescopes.
HST and JWST SIAF entries are supported via class inheritance. Frames, transformations, conventions and property/attribute names are as defined for JWST in JWST-STScI-001550.
Transformations between four coordinate systems (“frames”) are supported:
Detector (“det”) : units of pixels, according to detector read out axes orientation as defined by SIAF. This system generally differs from the JWST-DMS detector frame definition.
Science (“sci”) : units of pixels, corresponds to DMS coordinate system.
Ideal (“idl”) : units of arcseconds, usually a tangent plane projection with reference point at aperture reference location.
Telescope or V2/V3 (“tel”) : units of arcsecs, spherical coordinate system.
Examples
extract one aperture from a Siaf object:
ap = some_siaf['desired_aperture_name']
convert pixel coordinates to V2V3 / tel coords (takes pixel coords, returns arcsec):
ap.det_to_tel(1024, 512)
convert idl coords to sci pixels:
ap.idl_to_sci( 10, 3)
Frames can also be defined by strings:
ap.convert(1024, 512, from_frame='det', to_frame='tel')
Get aperture vertices/corners in tel frame:
ap.corners('tel')
Get the reference point defined in the SIAF:
ap.reference_point('tel')
plot coordinates in idl frame:
ap.plot('tel')
There exist methods for all of the possible
{tel,idl,sci,det}_to_{tel,idl,sci,det}
combinations.Set attributes required for PRD.
- closed_polygon_points(to_frame, rederive=True)[source]¶
Compute closed polygon points of aperture outline. Used for plotting and path generation.
- Parameters:
- to_framestr
Name of frame.
- rederivebool
Whether to rederive vertices from scratch (if True) or use idl values stored in SIAF.
- Returns:
- points_x, points_ytuple of numpy arrays
x and y coordinates of aperture vertices
- complement()[source]¶
Complement the attributes of an aperture.
- This method is useful when generating new apertures. The attributes that will be added are:
X[Y]IdlVert1[2,3,4]
X[Y]SciScale
- convert(x, y, from_frame, to_frame)[source]¶
Convert input coordinates from one frame to another frame.
- Parameters:
- xfloat
first coordinate
- yfloat
second coordinate
- from_framestr
Frame in which x,y are given
- to_framestr
Frame to transform into
- Returns:
- x’, y’tuple
Coordinates in to_frame
- corners(to_frame, rederive=True)[source]¶
Return coordinates of the aperture vertices in the specified frame.
The positions refer to the outside corners of the corner pixels, not the pixel centers.
- Parameters:
- to_framestr
Frame of returned coordinates
- rederivebool
If True, the parameters given in the aperture definition file are used to derive the vertices. Otherwise computations are made on the basis of the IdlVert attributes.
- Returns:
- x, ytuple of float arrays
coordinates of vertices in to_frame
- correct_for_dva(v2_arcsec, v3_arcsec, verbose=False)[source]¶
Apply differential velocity aberration correction to input arrays of V2/V3 coordinates.
Currently only implemented for HST apertures.
- Parameters:
- v2_arcsecfloat
v2 coordinates
- v3_arcsecfloat
v3 coordinates
- verbosebool
verbosity
- Returns:
- v2,v3tuple ot floats
Corrected coordinates
- det_to_sci(x_det, y_det, *args)[source]¶
Detector to science frame transformation, following Section 4.1 of JWST-STScI-001550.
- detector_transform(from_system, to_system, angle_deg=None, parity=None)[source]¶
Generate transformation model to transform between det <-> sci.
DetSciYAngle_deg, DetSciParity can be defined externally to override aperture attribute.
- Parameters:
- from_systemstr
Originating system
- to_systemstr
System to transform into
- angle_degfloat, int
DetSciYAngle
- parityint
DetSciParity
- Returns:
- x_model, y_modeltuple of astropy.modeling models
- distortion_transform(from_system, to_system, include_offset=True)[source]¶
Return polynomial distortion transformation between sci<->idl.
- Parameters:
- from_systemstr
Starting system (e.g. “sci”, “idl”)
- to_systemstr
Ending coordinate system (e.g. “sci” , “idl”)
- include_offsetbool
Whether to include the zero-point offset.
- Returns:
- x_modelastropy.modeling.Model
Correction in x
- y_modelastropy.modeling.Model
Correction in y
- dms_corner()[source]¶
Compute the pixel value of the lower left corner of an aperture.
This corresponds to the OSS corner position (x: ColCorner, y: RowCorner). The notation for OSS is 1-based, i.e. the lower left corner of a FULL subarray is (1,1)
- get_polynomial_coefficients()[source]¶
Return a dictionary of arrays holding the significant idl/sci coefficients.
- get_polynomial_derivatives(location='fiducial', coefficient_seed='Sci2Idl')[source]¶
Return derivative values for scale and rotation derivation.
The four partial derivatives of coefficients that transform between two frames E and R are:
b=(dx_E)/(dx_R ) c=(dx_E)/(dy_R ) e=(dy_E)/(dx_R ) f=(dy_E)/(dy_R )
- Parameters:
- locationstr or dict
if dict, has to have ‘x’ and ‘y’ elements
- coefficient_seed
- Returns:
- dict
- get_polynomial_linear_parameters(location='fiducial', coefficient_seed='Sci2Idl')[source]¶
Return linear polynomial parameters.
- Parameters:
- location
- coefficient_seed
- Returns:
- idl_to_tel(x_idl, y_idl, V3IdlYAngle_deg=None, V2Ref_arcsec=None, V3Ref_arcsec=None, method='planar_approximation', input_coordinates='tangent_plane', output_coordinates='tangent_plane', verbose=False)[source]¶
Convert from ideal to telescope (V2/V3) coordinate system.
By default, this implements the planar approximation, which is adequate for most purposes but may not be for all. Error is about 1.7 mas at 10 arcminutes from the tangent point. See JWST-STScI-1550 for more details. For higher accuracy, set method=’spherical_transformation’ in which case 3D matrix rotations are applied. Also by default, the input coordinates are in a tangent plane with a reference points at the origin (0,0) of the ideal frame. input_coordinates can be set to ‘spherical’.
- Parameters:
- x_idl: float
ideal X coordinate in arcsec
- y_idl: float
ideal Y coordinate in arcsec
- V3IdlYAngle_deg: float
overwrites self.V3IdlYAngle
- V2Ref_arcsecfloat
overwrites self.V2Ref
- V3Ref_arcsecfloat
overwrites self.V3Ref
- methodstr
must be one of [‘planar_approximation’, ‘spherical_transformation’, ‘spherical’]
- input_coordinatesstr
must be one of [‘tangent_plane’, ‘spherical’, ‘planar’]
- Returns:
- tuple of floats containing V2, V3 coordinates in arcsec
- path(to_frame)[source]¶
Return matplotlib path generated from aperture vertices.
- Parameters:
- to_framestr
Coordinate frame.
- Returns:
- pathmatplotlib.path.Path object
Path describing the aperture outline
- plot(frame='tel', label=False, ax=None, title=False, units='arcsec', show_frame_origin=None, mark_ref=False, fill=True, fill_color='cyan', fill_alpha=None, label_rotation=0.0, **kwargs)[source]¶
Plot this aperture.
- Parameters:
- framestr
Which coordinate system to plot in: ‘tel’, ‘idl’, ‘sci’, ‘det’
- labelstr or bool
Add text label. If True, text will be the default aperture name.
- axmatplotlib.Axes
Desired destination axes to plot into (If None, current axes are inferred)
- titlestr
Figure title. If set, add a label to the plot indicating which frame was plotted.
- unitsstr
one of ‘arcsec’, ‘arcmin’, ‘deg’
- show_frame_originstr or list
Plot frame origin (goes to plot_frame_origin()): None, ‘all’, ‘det’, ‘sci’, ‘raw’, ‘idl’, or a list of these.
- mark_refbool
Add marker for the (V2Ref, V3Ref) reference point defining this aperture.
- fillbool
Whether to color fill the aperture
- fill_colorstr
Fill color
- fill_alphafloat
alpha parameter for filled aperture
- colormatplotlib-compatible color
Color specification for this aperture’s outline, passed through to matplotlib.Axes.plot
- kwargsdict
Dictionary of optional parameters passed to matplotlib
- .. todo:: plotting in Sky frame, requires attitude
- #elif system == ‘radec’:
- # if attitude_ref is not None:
- # vertices = np.array(
- # siaf_rotations.pointing(attitude_ref, self.points_closed.T[0],
- #self.points_closed.T[1])).T
- # self.path_radec = matplotlib.path.Path(vertices)
- # else:
- # error(‘Please provide attitude_ref’)
- plot_detector_channels(frame, color='0.5', alpha=0.3, evenoddratio=0.5, ax=None)[source]¶
Outline the detector readout channels.
These are depicted as alternating light/dark bars to show the regions read out by each of the output amps.
- Parameters:
- framestr
Which coordinate system to plot in: “tel”, “idl”, “sci”, “det”
- colormatplotlib-compatible color
Color specification for the amplifier shaded region, passed through to matplotlib.patches.Polygon as facecolor
- alphafloat
Opacity of odd-numbered amplifier region overlays (for even, see “evenoddratio”)
- evenoddratiofloat
Ratio of opacity between even and odd amplifier region overlays
- axmatplotlib.Axes
Desired destination axes to plot into (If None, current axes are inferred from pyplot.)
- plot_frame_origin(frame, which='sci', units='arcsec', ax=None)[source]¶
Indicate the origins of the detector and science frames.
Red and blue squares indicate the detector and science frame origin, respectively.
- Parameters:
- framestr
Which coordinate system to plot in: ‘tel’, ‘idl’, ‘sci’, ‘det’
- whichstr of list
Which origin to plot: ‘all’, ‘det’, ‘sci’, ‘raw’, ‘idl’, or a list
- unitsstr
one of ‘arcsec’, ‘arcmin’, ‘deg’
- axmatplotlib.Axes
Desired destination axes to plot into (If None, current axes are inferred from pyplot.)
- raw_to_sci(x_raw, y_raw)[source]¶
Convert from raw/native coordinates to SIAF-Science coordinates.
SIAF-Science coordinates are the same as DMS coordinates for FULLSCA apertures. Implements the fits_generator description described the table attached to https://jira.stsci.edu/browse/JWSTSIAF-25 (except for Guider 2) and implemented in https://github.com/STScI-JWST/jwst/blob/master/jwst/fits_generator/create_dms_data.py
- Parameters:
- x_rawfloat
Raw x coordinate
- y_rawfloat
Raw y coordinate
- Returns:
- x_sci, y_scituple of science coordinates
References
See https://jwst-docs.stsci.edu/display/JDAT/Coordinate+Systems+and+Transformations See also JWST-STScI-002566 and JWST-STScI-003222 Rev A
- sci_to_det(x_sci, y_sci, *args)[source]¶
Science to detector frame transformation, following Section 4.1 of JWST-STScI-001550.
- sci_to_raw(x_sci, y_sci)[source]¶
Convert from Science coordinates to raw/native coordinates.
Implements the fits_generator description described the table attached to https://jira.stsci.edu/browse/JWSTSIAF-25 (except for Guider 2) and implemented in https://github.com/STScI-JWST/jwst/blob/master/jwst/fits_generator/create_dms_data.py
- Parameters:
- x_scifloat
Science x coordinate
- y_scifloat
Science y coordinate
- Returns:
- x_raw, y_rawtuple of raw coordinates
References
See https://jwst-docs.stsci.edu/display/JDAT/Coordinate+Systems+and+Transformations See also JWST-STScI-002566 and JWST-STScI-003222 Rev A
- set_attitude_matrix(attmat)[source]¶
Set an atttude matrix, for use in subsequent transforms to sky frame.
- set_distortion_coefficients_from_file(file_name)[source]¶
Set polynomial from standardized file.
- Parameters:
- file_namestr
Reference file containing coefficients
- Returns:
- set_polynomial_coefficients(sci2idlx, sci2idly, idl2scix, idl2sciy)[source]¶
Set the values of polynomial coefficients.
- Parameters:
- sci2idlxarray
Coefficients
- sci2idlyarray
Coefficients
- idl2scixarray
Coefficients
- idl2sciyarray
Coefficients
- sky_to_tel(*args)[source]¶
Sky to Tel frame transformation. Requires an attitude matrix.
Input sky coords should be given in degrees RA and Dec, or equivalent astropy.Quantities
Results are returned as floats with implicit units of arcseconds, for consistency with the other *_to_tel functions.
- tel_to_idl(v2_arcsec, v3_arcsec, V3IdlYAngle_deg=None, V2Ref_arcsec=None, V3Ref_arcsec=None, method='planar_approximation', output_coordinates='tangent_plane', input_coordinates='tangent_plane')[source]¶
Convert from telescope (V2/V3) to ideal coordinate system.
By default, this implements the planar approximation, which is adequate for most purposes but may not be for all. Error is about 1.7 mas at 10 arcminutes from the tangent point. See JWST-STScI-1550 for more details. For higher accuracy, set method=’spherical_transformation’ in which case 3D matrix rotations are applied.
Also by default, the output coordinates are in a tangent plane with a reference point at the origin (0,0) of the ideal frame.
- Parameters:
- v2_arcsecfloat
V2 coordinate in arcsec
- v3_arcsecfloat
V3 coordinate in arcsec
- V3IdlYAngle_degfloat
overwrites self.V3IdlYAngle
- V2Ref_arcsecfloat
overwrites self.V2Ref
- V3Ref_arcsecfloat
overwrites self.V3Ref
- methodstr
must be one of [‘planar_approximation’, ‘spherical_transformation’]
- output_coordinatesstr
must be one of [‘tangent_plane’, ‘spherical’]
- Returns:
- tuple of floats containing x_idl, y_idl coordinates in arcsec
- tel_to_sky(*args)[source]¶
Tel to sky frame transformation. Requires an attitude matrix.
Note, sky coordinates are returned as floats with angles in degrees. It’s easy to cast into an astropy SkyCoord and then convert into any desired other representation. For instance:
>>> sc = astropy.coordinates.SkyCoord( *aperture.tel_to_sky(xtel, ytel), unit='deg' ) >>> sc.to_string('hmsdms')
- telescope_transform(from_system, to_system, V3IdlYAngle_deg=None, V2Ref_arcsec=None, V3Ref_arcsec=None, verbose=False)[source]¶
Return transformation model between tel<->idl.
V3IdlYAngle_deg, V2Ref_arcsec, V3Ref_arcsec can be given as arguments to override aperture attributes.
- Parameters:
- from_systemstr
Starting system (“tel” or “idl”)
- to_systemstr
Ending coordinate system (“tel” or “idl”)
- V3IdlYAngle_degfloat
Angle
- V2Ref_arcsecfloat
Reference coordinate
- V3Ref_arcsecfloat
Reference coordinate
- verbosebool
verbosity
- Returns:
- x_model, y_modeltuple of astropy.modeling models
- class pysiaf.aperture.HstAperture[source]¶
Class for apertures of HST instruments.
Initialize HST aperture by inheriting from Aperture class.
- closed_polygon_points(to_frame, rederive=False)[source]¶
Compute closed polygon points of aperture outline.
- compute_tvs_matrix(v2_arcsec=None, v3_arcsec=None, pa_deg=None, verbose=False)[source]¶
Compute the TVS matrix from tvs-specific v2,v3,pa parameters.
- Parameters:
- v2_arcsecfloat
offset
- v3_arcsecfloat
offset
- pa_degfloat
angle
- verbosebool
verbosity
- Returns:
- tvsnumpy matrix
TVS matrix
- corners(to_frame, rederive=False)[source]¶
Return coordinates of the aperture vertices in the specified frame.
- idl_to_tel(x_idl, y_idl, V3IdlYAngle_deg=None, V2Ref_arcsec=None, V3Ref_arcsec=None, method='planar_approximation', input_coordinates='tangent_plane', output_coordinates=None, verbose=False)[source]¶
Convert ideal coordinates to telescope (V2/V3) coordinates for HST.
INPUT COORDINATES HAVE TO BE FGS OBJECT SPACE CARTESIAN X,Y coordinates
For HST FGS, transformation is implemented using the FGS TVS matrix. Parameter names are overloaded for simplicity: tvs_pa_deg = V3IdlYAngle_deg tvs_v2_arcsec = V2Ref_arcsec tvs_v3_arcsec = V3Ref_arcsec
- Parameters:
- x_idlfloat
ideal x coordinates in arcsec
- y_idlfloat
ideal y coordinates in arcsec
- V3IdlYAngle_degfloat
angle
- V2Ref_arcsecfloat
Reference coordinate
- V3Ref_arcsecfloat
Reference coordinate
- methodstr
- input_coordinatesstr
- Returns:
- v2, v3tuple of float
Telescope coordinates in arcsec
- rearrange_fgs_alignment_parameters(pa_deg_in, v2_arcsec_in, v3_arcsec_in, direction)[source]¶
Convert to/from alignment parameters make FGS look like a regular camera aperture.
- Parameters:
- pa_deg_infloat
- v2_arcsec_infloat
- v3_arcsec_infloat
- directionstr
One of [‘fgs_to_camera’, ]
- Returns:
- pa, v2, v3tuple
re-arranged and sometimes sign-inverted alignment parameters
- set_idl_reference_point(v2_ref, v3_ref, verbose=False)[source]¶
Determine the reference point in the Ideal frame.
V2Ref and V3Ref are determined via the TVS matrix. The tvs parameters that determine the TVS matrix itself are derived and added to the attribute list.
- Parameters:
- v2_reffloat
reference coordinate
- v3_reffloat
reference coordinate
- verbosebool
verbosity
- set_tel_reference_point(verbose=True)[source]¶
Recompute and set V2Ref and V3Ref to actual position in tel/V frame.
This is after using those V2Ref and V3Ref attributes for TVS matrix update.
- tel_to_idl(v2_arcsec, v3_arcsec, V3IdlYAngle_deg=None, V2Ref_arcsec=None, V3Ref_arcsec=None, method='planar_approximation', output_coordinates='tangent_plane', input_coordinates='tangent_plane')[source]¶
Convert from telescope (V2/V3) to ideal coordinate system.
By default, this implements the planar approximation, which is adequate for most purposes but may not be for all. Error is about 1.7 mas at 10 arcminutes from the tangent point. See JWST-STScI-1550 for more details. For higher accuracy, set method=’spherical_transformation’ in which case 3D matrix rotations are applied.
Also by default, the output coordinates are in a tangent plane with a reference point at the origin (0,0) of the ideal frame.
- Parameters:
- v2_arcsecfloat
V2 coordinate in arcsec
- v3_arcsecfloat
V3 coordinate in arcsec
- V3IdlYAngle_degfloat
overwrites self.V3IdlYAngle
- V2Ref_arcsecfloat
overwrites self.V2Ref
- V3Ref_arcsecfloat
overwrites self.V3Ref
- methodstr
must be one of [‘planar_approximation’, ‘spherical_transformation’]
- output_coordinatesstr
must be one of [‘tangent_plane’, ‘spherical’]
- Returns:
- tuple of floats containing x_idl, y_idl coordinates in arcsec
- class pysiaf.aperture.JwstAperture[source]¶
Class for apertures of JWST instruments.
Initialize JWST aperture by inheriting from Aperture.
- class pysiaf.aperture.NirspecAperture(tilt=None, filter_name='CLEAR')[source]¶
Class for apertures of the JWST NIRSpec instrument.
Initialize NIRSpec aperture by inheriting from JWSTAperture.
- det_to_sci(x_det, y_det, *args)[source]¶
Detector to science frame transformation. Use parent aperture if SLIT.
- gwa_to_ote(gwa_x, gwa_y)[source]¶
NIRSpec transformation from GWA sky side to OTE frame XAN, YAN.
- Parameters:
- gwa_xfloat
GWA coordinate
- gwa_yfloat
GWA coordinate
- Returns:
- xan, yantuple of floats
XAN, YAN in degrees
- gwa_to_sci(x_gwa, y_gwa)[source]¶
NIRSpec transformation from GWA detector side to Science frame.
- Parameters:
- x_gwafloat
GWA coordinate
- y_scifloat
GWA coordinate
- Returns:
- x_sci, y_scituple of floats
Science coordinates
- gwain_to_gwaout(x_gwa, y_gwa)[source]¶
Transform from GWA detector side to GWA skyward side. This is the effect of the mirror.
- Parameters:
- x_gwafloat
GWA coordinate
- y_gwafloat
GWA coordinate
- Returns:
- x_gwap, y_gwaptuple of floats
GWA skyward side coordinates
References
Giardino, Ferruit, and Alves de Oliveira (2014), ESA NTN-2014-005 Proffitt et al. (2018), JWST-STScI-005921
- gwaout_to_gwain(x_gwa, y_gwa)[source]¶
Transform from GWA skyward side to GWA detector side. Effect of mirror.
- Parameters:
- x_gwafloat
GWA coordinate
- y_gwafloat
GWA coordinate
- Returns:
- x_gwap, y_gwaptuple of floats
GWA detector side coordinates
References
Equations for the reverse transform T. Keyes (private communication) will be documented in next update of JWST-STScI-005921.
- idl_to_sci(x_idl, y_idl)[source]¶
Transform to Science frame using special implementation for NIRSpec via tel frame.
- Parameters:
- x_idlfloat
Ideal coordinate
- y_idlfloat
Ideal coordinate
- Returns:
- x_sci, y_scituple of floats
Science coordinates
- ote_to_gwa(ote_x, ote_y)[source]¶
NIRSpec transformation from OTE frame XAN, YAN to GWA sky side.
- Parameters:
- ote_xfloat
XAN coordinate in degrees
- ote_yfloat
YAN coordinate in degrees
- Returns:
- gwa_x, gwa_ytuple of floats
GWA coordinates
- sci_to_det(x_sci, y_sci, *args)[source]¶
Science to detector frame transformation. Use parent aperture if SLIT.
- sci_to_gwa(x_sci, y_sci)[source]¶
NIRSpec transformation from Science frame to GWA detector side.
- Parameters:
- x_scifloat
Science frame coordinate
- y_scifloat
Science frame coordinate
- Returns:
- x_gwa, y_gwatuple of floats
GWA detector side coordinates
- sci_to_idl(x_sci, y_sci)[source]¶
Transform to ideal frame, using special implementation for NIRSpec via tel frame.
- Parameters:
- x_scifloat
Science coordinate
- y_scifloat
Science coordinate
- Returns:
- v2, v3tuple of floats
telescope coordinates
- class pysiaf.aperture.RomanAperture[source]¶
A class that defines the accepted aperture types for the Roman SIAF. This is built upon the pysiaf.aperture.Aperture base class.
Set attributes required for PRD.
- class pysiaf.aperture.SiafCoordinates(x, y, frame)[source]¶
A helper class to hold coordinate arrays associated with a SIAF frame.
Initialize object.
- pysiaf.aperture.compare_apertures(reference_aperture, comparison_aperture, absolute_tolerance=None, attribute_list=None, print_file=<_io.TextIOWrapper name='<stdout>' mode='w' encoding='utf-8'>, fractional_tolerance=1e-06, verbose=False, ignore_attributes=None)[source]¶
Compare the attributes of two apertures.
- Parameters:
- reference_aperture
- comparison_aperture
- absolute_tolerance
- attribute_list
- print_file
- fractional_tolerance
- verbose
- ignore_attributes
- pysiaf.aperture.get_hst_to_jwst_coefficient_order(polynomial_degree)[source]¶
Return array of indices that convert an aeeay of HST coefficients to JWST ordering.
This assumes that the coefficient orders are as follows HST: 1, y, x, y^2, xy, x^2, y^3, xy^2, x^2y, x^3 … (according to Cox’s /grp/hst/OTA/alignment/FocalPlane.xlsx) JWST: 1, x, y, x^2, xy, y^2, x^3, x^2y, xy^2, y^3 …
- Parameters:
- polynomial_degreeint
polynomial degree
- Returns:
- conversion_indexnumpy array
Array of indices
- pysiaf.aperture.linear_transform_model(from_system, to_system, parity, angle_deg)[source]¶
Create an astropy.modeling.Model object for linear transforms: det <-> sci and idl <-> tel.
See sics_to_v2v3 (HST): x_v2v3 = v2_origin + parity * x_sics * np.cos(np.deg2rad(theta_deg)) + y_sics * np.sin( np.deg2rad(theta_deg)) y_v2v3 = v3_origin - parity * x_sics * np.sin(np.deg2rad(theta_deg)) + y_sics * np.cos( np.deg2rad(theta_deg))
- Parameters:
- from_systemstr
Starting system
- to_systemstr
End system
- parityint
Parity
- angle_degfloat
Angle in degrees
- Returns:
- xmodel, ymodeltuple of astropy models
Transformation models
- pysiaf.aperture.points_on_arc(x0, y0, radius, phi1_deg, phi2_deg, N=100)[source]¶
Compute points that lie on a circular arc between angles phi1 and phi2.
- Parameters:
- x0float
offset
- y0float
offset
- radiusfloat
Radius of circle
- phi1_degfloat
Start angle of pie in deg
- phi2_degfloat
End angle of pie in deg
- Nint
Number of sampled points
- Returns:
- x, ytuple of floats
Coordinates on circle
- pysiaf.aperture.to_distortion_model(coefficients, degree=5)[source]¶
Create an astropy.modeling.Model object for distortion polynomial.
- Parameters:
- coefficientsarray like
Coefficients from the ISIM transformations file.
- degreeint
Degree of polynomial. Default is 5 as in the ISIM file but many of the polynomials are of a smaller degree.
- Returns:
- polyastropy.modeling.Polynomial2D
Polynomial model transforming one coordinate (x or y) between two systems.
compare¶
Functions to support comparisons between SIAF files.
References¶
The format of the difference files was adapted from Colin Cox’ matchcsv.py dict_compare was adapted from https://stackoverflow.com/questions/4527942/comparing-two-dictionaries-in-python/4527957
- pysiaf.utils.compare.compare_inspection_figures(comparison_siaf_input, reference_siaf_input=None, report_dir=None, selected_aperture_name=None, skipped_aperture_type=None, tags=None, mark_ref=False, xlimits=None, ylimits=None, filename_appendix='', label=False)[source]¶
Visualize aperture of two SIAF files.
- Parameters:
- comparison_siaf_inputstr (absolute file name) or pysiaf.Siaf object
The SIAF that will be compared to the reference_siaf
- reference_siaf_inputstr (absolute file name) or pysiaf.Siaf object
The reference SIAF. Defaults to the current PRD content.
- report_dir
- selected_aperture_name str or list
aperture name(s) to include in plot
- skipped_aperture_typestr or list
aperture type(s) not to include in plot
- tags
- xlimitstuple of limits of output plots
- ylimitstuple of limits of output plots
- Returns:
- pysiaf.utils.compare.compare_siaf(comparison_siaf_input, fractional_tolerance=0.0001, reference_siaf_input=None, report_file=None, report_dir=None, verbose=True, make_figures=False, selected_aperture_name=None, ignore_attributes=None, tags=None)[source]¶
Compare two SIAF files and write a difference file. Generate comparison figures showing the apertures if specified.
- Parameters:
- comparison_siaf_inputstr (absolute file name) or pysiaf.Siaf object
The SIAF that will be compared to the reference_siaf
- fractional_tolerancefloat
Value above which a fractional difference in parameter value will be written in the report
- reference_siaf_inputstr (absolute file name) or pysiaf.Siaf object
The reference SIAF. Defaults to the current PRD content.
- report_file
- report_dir
- verbose
- make_figures
- selected_aperture_name
- pysiaf.utils.compare.compare_transformation_roundtrip(comparison_siaf_input, fractional_tolerance=0.0001, reference_siaf_input=None, report_file=None, report_dir=None, verbose=True, make_figures=False, selected_aperture_name=None, skipped_aperture_type=None, instrument=None, make_plot=False, tags=None)[source]¶
Compare the forward-backward roundtrip transformations of two SIAF files. and write a difference file.
- Parameters:
- comparison_siaf_inputstr (absolute file name) or pysiaf.Siaf object
The SIAF that will be compared to the reference_siaf
- fractional_tolerancefloat
Value above which a fractional difference in parameter value will be written in the report
- reference_siaf_inputstr (absolute file name) or pysiaf.Siaf object
The reference SIAF. Defaults to the current PRD content.
- report_file
- report_dir
- verbose
- make_figures
- selected_aperture_namestr or list
aperture name(s) to include in plot
- skipped_aperture_typestr or list
aperture type(s) not to include in plot
- instrument
- make_plot
- Returns:
- roundtrip_tableastropy.table.Table object
table containing roundtrip data
- pysiaf.utils.compare.dict_compare(dictionary_1, dictionary_2)[source]¶
Compare two dictionaries and return keys of the differing items. From https://stackoverflow.com/questions/4527942/comparing-two-dictionaries-in-python/4527957
- Parameters:
- dictionary_1dict
first dictionary
- dictionary_2dict
second dictionary
- Returns:
- added, removed, modified, sameset
Sets of dictionary keys that were added, removed, modified, or are the same
polynomial¶
A collection of functions to manipulate polynomials and their coefficients.
- pysiaf.utils.polynomial.add_rotation(A, B, theta_deg)[source]¶
Add rotation after polynomial transformation.
Use when a distortion transformation using polynomials A and B is followed by a rotation.
u = A(x,y) v = B(x,y) followed by u2 = u*cos(theta) + v*sin(theta) v2 = -u*sin(theta) + v*cos(theta) This routine supplies a modified pair of polynomial which combine both steps i.e u2 = A2(x,y), v2 = B2(x,y)
- Parameters:
- Aarray
Set of polynomial coefficients converting from (x,y) to a variable u
- Barray
set of polynomial coefficients converting from(x,y) to avariable v
- theta_degfloat
The angle in degrees of a rotationin the (u,v) plane
- Returns:
- A2array
set of polynomial coefficiients providing combined steps from (x,y) to u2
- B2array
set of polynomial coefficients providing combined steps from (x,y) to v2
Notes
Function formerly named Rotate or rotate_coefficients. Ported from makeSIAF.py by J. Sahlmann 2018-01-03.
- pysiaf.utils.polynomial.choose(n, r)[source]¶
Return number of ways of choosing r items from an array with n items.
- Parameters:
- nint
number of items to choose from
- rint
number of items to choose
- Returns:
- combinationsint
The number if ways of making the choice
- pysiaf.utils.polynomial.dpdx(a, x, y)[source]¶
Differential with respect to x.
The polynomial is defined as p(x,y) = a[i,j] * x**(i-j) * y**j summed over i and j. Then dp/dx = (i-j) * a[i,j] * x**(i-j-1) * y**j.
- Parameters:
- aarray
a linear array of polynomial coefficients in JWST order.
- xarray
an integer or float variable(or an array of same) representing pixel x positions
- yarray
a variable (or an array) representing pixel y positions. x and y must be of same shape.
- Returns:
- differentialarray
float values of dp/dx for the given (x,y) point(s)
- pysiaf.utils.polynomial.dpdy(a, x, y)[source]¶
Differential with respect to y.
The polynomial is defined as p(x,y) = a[i,j] * x**(i-j) * y**j, summed over i and j Then dp/dy = j * a[i,j] * x**(i-j) * y**(j-1)
- Parameters:
- aarray
an array of polynomial coefficients in JWST arrangement. The number of coefficients must be (order+1)(order+2)/2
- xarray
an integer or float variable(or an array of same) representing pixel x positions
- yarray
a variable (or an array) representing pixel y positions
- Returns:
- differentialarray
float value of dp/dy for the given (x,y) point(s) where p(x,y) is the value of the polynomial
- pysiaf.utils.polynomial.flatten(coefficients)[source]¶
Convert triangular layout to linear array.
For many of the polynomial operations the coefficients A[i,j] are contained in an array of dimension (order+1, order+1) but with all elements where j > i set equal to zero. This is what is called the triangular layout.
The flattened layout is a one-dimensional array containing only the elements where j <= i.
- Parameters:
- coefficientsarray
a two-dimensional float array of shape (order+1, order+1) supplying the polynomial coefficients in JWST order. The coefficient at position [i,j] multiplies the value of x**(i-j) * y**j.
- Returns:
- flate_coefficientsarray
a one-dimensional array including only those terms where i <= j
- pysiaf.utils.polynomial.flip_x(A)[source]¶
Change sign of all coefficients with odd x power.
Used when we have a polynomial expansion in terms of variables x and y and we wish to obtain one in which the sign of x is reversed
- Parameters:
- Aarray
A set of polynomial coefficients given in the triangular layout as described in poly
- Returns:
- AFarray
Modified or flipped set of coefficients matching negated x values.
- pysiaf.utils.polynomial.flip_xy(A)[source]¶
Change sign for coeffs where sum of x and y powers is odd.
Used when we have a polynomial expansion in terms of variables x and y and we wish to obtain one in which the signs of x and y are reversed
- Parameters:
- Aarray
A set of polynomial coefficients given in the triangular layout as described in the function poly
- Returns:
- AFarray
Modified or flipped set of coefficients matching negated x and y values.
- pysiaf.utils.polynomial.flip_y(A)[source]¶
Change sign of all coefficients with odd y power.
Used when we have a polynomial expansion in terms of variables x and y and we wish to obtain one in which the sign of y is reversed
- Parameters:
- Aarray
A set of polynomial coefficients given in the triangular layout as described in the function poly
- Returns:
- AFarray
Modified or flipped set of coefficients matching negated y values.
- pysiaf.utils.polynomial.invert(A, B, u, v, verbose=False)[source]¶
Newton Raphson method in two dimensions.
Given that u = A[i,j] * x**(i-j) * y**j and v = B[i,j] * x**(i-j) * y**j find the values of x and y from the values of u and v
- Parameters:
- Aarray
A set of polynomial coefficients given in the linear layout as described in the function poly converting (x,y) to u
- Barray
A set of polynomial coefficients given in the linear layout as described in the function poly converting (x,y) to v
- uarray
The result of applying the A coefficients to the (x,y) position
- varray
The result of applying the B coefficients to the (x, y)position
- verbosebool
Logical variable, set True if full text output required
- Returns:
- x, ytuple of arrays
The pair of values which transform to (u,v)
- errfloat
the standard deviation of the fit
- iterationint
the number of iterations taken to determine the solution
- pysiaf.utils.polynomial.jacob(a, b, x, y)[source]¶
Calculate relative area using the Jacobian.
` `
` | da_dx db_dx | `
`Jacobian = | | `
` | da_dy db_dy | `
` `
Then the relative area is the absolute value of the determinant of the Jacobian. x and y will usually be Science coordinates while u and v are Ideal coordinates
- Parameters:
- aarray
set of polynomial coefficients converting from (x,y) to u
- barray
set of polynomial coefficients converting from (x,y) to v
- xarray
x pixel position or array of x positions
- yarray
y pixel position or array of y positions matching the y positions
- Returns:
- areaarray
area in (u,v) coordinates matching unit area in the (x,y) coordinates.
- pysiaf.utils.polynomial.number_of_coefficients(poly_degree)[source]¶
Return number of coefficients corresponding to polynomial degree.
- pysiaf.utils.polynomial.poly(a, x, y, order=4)[source]¶
Polynomial evaluation.
pol = a[i,j] * x**(i-j) * y**j summed over i and j, where i runs from 0 to order. Then for each value of i, j runs from 0 to i. For many of the polynomial operations the coefficients A[i,j] are contained in an array of dimension (order+1, order+1) but with all elements where j > i set equal to zero. This is called the triangular layout. The flattened layout is a one-dimensional array containing copies of only the elements where j <= i.
The JWST layout is a[0,0] a[1,0] a[1,1] a[2,0] a[2,1] a[2,2] … The number of coefficients will be (n+1)(n+2)/2
- Parameters:
- aarray
float array of polynomial coefficients in flattened arrangement
- xarray
x pixel position. Can be integer or float or an array of integers or floats
- yarray
y pixel position in same layout as x positions.
- orderint
integer polynomial order
- Returns:
- polfloat
result as described above
- pysiaf.utils.polynomial.polyfit(u, x, y, order)[source]¶
Fit polynomial to a set of u values on an x,y grid.
u is a function u(x,y) being a polynomial of the form u = a[i, j] x**(i-j) y**j. x and y can be on a grid or be arbitrary values This version uses scipy.linalg.solve instead of matrix inversion. u, x and y must have the same shape and may be 2D grids of values.
- Parameters:
- uarray
an array of values to be the results of applying the sought after polynomial to the values (x,y)
- xarray
an array of x values
- yarray
an array of y values
- orderint
the polynomial order
- Returns:
- coeffs: array
polynomial coefficients being the solution to the fit.
- pysiaf.utils.polynomial.polynomial_degree(number_of_coefficients)[source]¶
Return degree of the polynomial that has number_of_coefficients.
- Parameters:
- number_of_coefficientsint
Number of polynomial coefficients
- Returns:
- polynomial_degreeint
Degree of the polynomial
- pysiaf.utils.polynomial.prepend_rotation_to_polynomial(a, theta, verbose=False)[source]¶
Rotate axes of coefficients by theta degrees.
Used when a distortion transformation using polynomials A and B is preceded by a rotation. The set of polynomial coefficients a[i,j] transform (x,y) as u = a[i,j] * x**(i-j) * y**j Summation over repeated indices is implied. If now we have a set of variables (xp,yp) rotated from (x,y) so that xp = x * cos(theta) - y * sin(theta) yp = x * sin(theta) + y * cos(theta) find a set of polynomial coefficients ap so that the same value of u is obtained from (xp,yp) i.e, u = ap[i,j] * xp**(i-j) * yp**j The rotation is opposite to the usual rotation as this routine was designed for the inverse transformation between Ideal and V2V3 or tel. Effectively the angle is reversed
- Parameters:
- aarray
Set of polynomial coefficients
- thetafloat
rotation angle in degrees
- verbosebool
logical variable set True only if print-out of coefficient factors is desired.
- Returns:
- arotatearray
set of coefficients derived as described above.
Notes
Function was formerly named RotateCoeffs.
- pysiaf.utils.polynomial.print_triangle(coefficients)[source]¶
Print coefficients in triangular layout.
A[0] A[1] A[2] A[3] A[4] A[5] … equivalent to A[0,0] A[1,0] A[1,1] A[2,0] A[2,1] A[2,2] … in [i,j] terms.
See method poly for details. This is just to display the coefficients. No calculation performed.
- Parameters:
- coefficientsarray
polynomial float array in linear layout
- pysiaf.utils.polynomial.reorder(A, B)[source]¶
Change coefficient order from y**2 xy x**2 to x**2 xy y**2 in both A and B.
- Parameters:
- Aarray
polynomial coefficients
- Barray
polynomial coefficients
- Returns:
- A2, B2: numpy arrays
coefficients with changed order
- pysiaf.utils.polynomial.rescale(A, B, C, D, scale)[source]¶
Apply a scale to forward and inverse polynomial coefficients.
- Parameters:
- Aarray
Polynomial coefficients
- Barray
Polynomial coefficients
- Carray
Polynomial coefficients
- Darray
Polynomial coefficients
- scalefloat
Scale factor to apply
- Returns:
- A_scaled, B_scaled, C_scaled, D_scaledtuple of numpy arrays
the scales coefficients
Notes
Ported from makeSIAF.py by J. Sahlmann 2018-01-03. J. Sahlmann 2018-01-04: fixed side-effect on ABCD variables
- pysiaf.utils.polynomial.rotation_scale_skew_from_derivatives(b, c, e, f)[source]¶
Compute rotations, scales, and skews from polynomial derivatives.
The four partial derivatives of coefficients that transform between two frames E and R are:
- Parameters:
- bfloat
(dx_E)/(dx_R )
- cfloat
(dx_E)/(dy_R )
- efloat
(dy_E)/(dx_R )
- ffloat
(dy_E)/(dy_R )
- Returns:
- pysiaf.utils.polynomial.shift_coefficients(a, xshift, yshift, verbose=False)[source]¶
Calculate coefficients of polynomial when shifted to new origin.
Given a polynomial function such that u = a[i,j] * x**[i-j] * y**[j] summed over i and j. Find the polynomial function ashift centered at xshift, yshift i.e the same value of u = ashift[i,j] * (x-xshift)**(i-j) * (y-yshift)**j.
- Parameters:
- aarray
Set of coefficients for a polynomial of the given order in JWST order
- xshiftfloat
x position in pixels of new solution center
- yshiftfloat
y position in pixels of new solution center
- verbosebool
logical variable to choose print-out of coefficient table - defaults to False
- Returns:
- ashiftarray
shifted version of the polynomial coefficients.
- pysiaf.utils.polynomial.transform_coefficients(A, a, b, c, d, verbose=False)[source]¶
Transform polynomial coefficients.
This allows for xp = a*x + b*y yp = c*x + d*y
- Parameters:
- Aarray
Polynomial coefficients
- afloat
factor
- bfloat
factor
- cfloat
factor
- dfloat
factor
- verbosebool
verbosity
- Returns:
- ATarray
Transformed coefficients
Notes
Designed to work with Sabatke solutions which included a linear transformation of the pixel coordinates before the polynomial distortion solution was calculated.
transform_coefficients
combines the two steps into a single polynomial.
- pysiaf.utils.polynomial.triangular_layout(coefficients)[source]¶
Convert linear array to 2-D array with triangular coefficient layout.
This is the reverse of the flatten method.
- Parameters:
- coefficientsarray
float array of polynomial coefficients. Must be of dimension (order+1)(order+2)/2.
- Returns:
- triangular_coefficientsarray
coefficients in triangular layout as described in poly.
- pysiaf.utils.polynomial.two_step(A, B, a, b)[source]¶
Combine linear step followed by a polynomial step into a single polynomial.
Designed to process Sabatke polynomials which had a linear transformation ahead of the polynomial fits.
Starting from a pair of polynomial arrays A and B such that u = A[i,j[ * xp**(i-j) * yp**j v = B[i,j] * xp**(i-j) * yp**j in which xp = a[0] + a[1].x + a[2].y yp = b[0] + b[1].x + b[2].y find AP and BP such that the same u and v values are given by u = AP[i,j] * x**(i-j) * y**j v = BP[i,j] * x**(i-j) * y**j All input and output polynomials are flattened arrays of dimension (order+1)(order+2)/2 Internally they are processed as equivalent two dimensional arrays as described in the poly documentation.
- Parameters:
- Aarray
polynomial converting from secondary xp and yp pixel positions to final coordinates u
- Barray
polynomial converting from secondary xp and yp pixel positions to final coordinates v
- Aflatarray
set of linear coefficients converting (x,y) to xp
- Bflatarray
set of linear coefficients converting (x,y) to yp
- Returns:
- Aflat, Bflattuple of arrays
polynomial coefficients as calculated
projection¶
A collection of functions to support tangent-plane de-/projections.
- pysiaf.utils.projection.deproject_from_tangent_plane(x, y, ra_ref, dec_ref, scale=1.0, unwrap=True)[source]¶
Convert pixel coordinates into ra/dec coordinates using a tangent plane de-projection.
The projection’s reference point has to be specified. See the inverse transformation radec2Pix_TAN.
- Parameters:
- xfloat or array of floats
Pixel coordinate (default is in decimal degrees, but depends on value of scale parameter) x/scale has to be degrees.
- yfloat or array of floats
Pixel coordinate (default is in decimal degrees, but depends on value of scale parameter) x/scale has to be degrees.
- ra_reffloat
Right Ascension of reference point in decimal degrees
- dec_ref: float
declination of reference point in decimal degrees
- scalefloat
Multiplicative factor that is applied to the input values. Default is 1.0
- Returns:
- rafloat
Right Ascension in decimal degrees
- dec: float
declination in decimal degrees
- pysiaf.utils.projection.project_to_tangent_plane(ra, dec, ra_ref, dec_ref, scale=1.0)[source]¶
Convert ra/dec coordinates into pixel coordinates using a tangent plane projection.
Theprojection’s reference point has to be specified. Scale is a convenience parameter that defaults to 1.0, in which case the returned pixel coordinates are also in degree. Scale can be set to a pixel scale to return detector coordinates in pixels
- Parameters:
- rafloat
Right Ascension in decimal degrees
- dec: float
declination in decimal degrees
- ra_reffloat
Right Ascension of reference point in decimal degrees
- dec_ref: float
declination of reference point in decimal degrees
- scalefloat
Multiplicative factor that is applied to the returned values. Default is 1.0
- Returns:
- x,yfloat
pixel coordinates in decimal degrees if scale = 1.0
read¶
Functions to read Science Instrument Aperture Files (SIAF) and SIAF reference files.
For JWST SIAF, reading XML and CSV format are supported. For HST SIAF, only .dat files can be read.
References¶
Parts of read_hst_siaf were adapted from Matt Lallo’s plotap.f. Parts of read_jwst_siaf were adapted from jwxml.
- pysiaf.iando.read.get_jwst_siaf_instrument(tree)[source]¶
Return the instrument specified in the first aperture of a SIAF xml tree.
- Returns:
- instrumentstr
All Caps instrument name, e.g. NIRSPEC
- pysiaf.iando.read.get_siaf(input_siaf, observatory='JWST')[source]¶
Return a Siaf object corresponding to input_siaf which can be a string path or a Siaf object.
- Parameters:
- input_siaf
- observatory
- Returns:
- siaf_object: pysiaf.Siaf
Siaf object
- pysiaf.iando.read.read_hst_fgs_amudotrep(file=None, version=None)[source]¶
Read HST FGS amu.rep file which contain the TVS matrices.
- Parameters:
- filepathstr
Path to file.
- Returns:
- datadict
Dictionary that holds the file content ordered by FGS number
- pysiaf.iando.read.read_hst_siaf(file=None, version=None)[source]¶
Read apertures from HST SIAF file and return a collection.
This was partially ported from Lallo’s plotap.f.
- Parameters:
- filestr
- AperNamesstr list
- Returns:
- apertures: dict
Dictionary of apertures
- pysiaf.iando.read.read_jwst_siaf(instrument=None, filename=None, basepath=None)[source]¶
Read the JWST SIAF and return a collection of apertures.
- Parameters:
- instrumentstr
instrument name (case-insensitive)
- filenamestr
Absolute path to alternative SIAF xml file
- basepathstr
Directory containing alternative SIAF xml file conforming with standard naming convention
- Returns:
- aperturesdict
dictionary of apertures
- pysiaf.iando.read.read_roman_siaf(siaf_file=None)[source]¶
- Returns:
- apertures (collections.OrderedDict)
An ordered dictionary of pysiaf.aperture objects containing each aperture from the Roman SIAF file that was read.
- pysiaf.iando.read.read_siaf_alignment_parameters(instrument)[source]¶
Return astropy table.
- Parameters:
- instrument
- Returns:
- : astropy table
- pysiaf.iando.read.read_siaf_aperture_definitions(instrument, directory=None)[source]¶
Return astropy table.
- Parameters:
- instrumentstr
instrument name (case insensitive)
- Returns:
- : astropy table
content of SIAF reference file
- pysiaf.iando.read.read_siaf_ddc_mapping_reference_file(instrument)[source]¶
Return dictionary with the DDC mapping.
- Parameters:
- instrumentstr
instrument name (case insensitive)
- Returns:
- : astropy table
- pysiaf.iando.read.read_siaf_detector_layout()[source]¶
Return the SIAF detector layout read from the SIAF reference file.
- Returns:
- : astropy table
- pysiaf.iando.read.read_siaf_detector_reference_file(instrument)[source]¶
Return astropy table.
- Parameters:
- instrumentstr
instrument name (case insensitive)
- Returns:
- : astropy table
- pysiaf.iando.read.read_siaf_distortion_coefficients(instrument=None, aperture_name=None, file_name=None)[source]¶
Return astropy table.
- Parameters:
- instrumentstr
instrument name (case insensitive)
- aperture_namestr
name of master aperture
- file_namestr
file name to read from, ignoring the first two arguments
- Returns:
- : astropy table
rotations¶
A collection of basic routines for performing rotation calculations.
- pysiaf.utils.rotations.attitude(v2, v3, ra, dec, pa)[source]¶
Return rotation matrix that transforms from v2,v3 to RA,Dec.
Makes a 3D rotation matrix which rotates a unit vector representing a v2,v3 position to a unit vector representing an RA, Dec pointing with an assigned position angle Described in JWST-STScI-001550, SM-12, section 5.1.
- Parameters:
- v2float
a position measured in arc-seconds
- v3float
a position measured in arc-seconds
- rafloat
Right Ascension on the sky in degrees
- decfloat
Declination on the sky in degrees
- pafloat
Position angle in degrees measured from North to V3 axis in North to East direction.
- Returns:
- mnumpy matrix
A (3 x 3) matrix represents the attitude of the telescope which points the given V2V3 position to the indicated RA and Dec and with the V3 axis rotated by position angle pa
- pysiaf.utils.rotations.attitude_matrix(nu2, nu3, ra, dec, pa, convention='JWST')[source]¶
Return attitude matrix.
Makes a 3D rotation matrix that transforms between telescope frame and sky. It rotates a unit vector on the idealized focal sphere (specified by the spherical coordinates nu2, nu3) to a unit vector representing an RA, Dec pointing with an assigned position angle measured at nu2, nu3. See JWST-STScI-001550, SM-12, section 5.1.
- Parameters:
- nu2float
an euler angle (default unit is arc-seconds)
- nu3float
an euler angle (default unit is arc-seconds)
- rafloat
Right Ascension on the sky in degrees
- decfloat
Declination on the sky in degrees
- pafloat
Position angle of V3 axis at nu2,nu3 measured from North to East (default unit is degree)
- Returns:
- mnumpy matrix
the attitude matrix
- pysiaf.utils.rotations.axial_rotation(ax, phi, vector)[source]¶
Apply direct rotation to a vector using Rodrigues’ formula.
- Parameters:
- axfloat array of size 3
a unit vector represent a rotation axis
- phifloat
angle in degrees to rotate original vector
- vectorfloat
array of size 3 representing any vector
- Returns:
- vfloat
array of size 3 representing the rotated vectot
- pysiaf.utils.rotations.convert_quantity(x_in, to_unit, factor=1.0)[source]¶
Check if astropy quantity and apply conversion factor
- Parameters:
- x_infloat or quantity
input
- to_unitastropy.units unit
unit to convert to
- factorfloat
Factor to apply if input is not a quantity
- Returns:
- x_outfloat
converted value
- pysiaf.utils.rotations.cross(a, b)[source]¶
Return cross product of two vectors c = a X b.
The order is significant. Reversing the order changes the sign of the result.
- Parameters:
- afloat array or list of length 3
first vector
- bfloat array or list of length 3
second vector
- Returns:
- c float array of length 3
the product vector
- pysiaf.utils.rotations.getv2v3(attitude, ra, dec)[source]¶
Return v2,v3 position of any RA and Dec using the inverse of attitude matrix.
- Parameters:
- attitude3 by 3 float array
the telescope attitude matrix
- rafloat
RA of sky position
- decfloat
Dec of sky position
- Returns:
- v2,v3tuple of floats
V2,V3 value at matching position
- pysiaf.utils.rotations.idl_to_tel_rotation_matrix(V2Ref_arcsec, V3Ref_arcsec, V3IdlYAngle_deg)[source]¶
Return 3D rotation matrix for ideal to telescope transformation.
- Parameters:
- V2Ref_arcsec
- V3Ref_arcsec
- V3IdlYAngle_deg
- Returns:
- pysiaf.utils.rotations.pointing(attitude, v2, v3, positive_ra=True, input_cartesian=False)[source]¶
Calculate where a v2v3 position points on the sky using the attitude matrix.
- Parameters:
- attitude3 by 3 float array
the telescope attitude matrix
- v2float or array of floats
V2 coordinate in arc-seconds
- v3float or array of floats
V3 coordinate in arc-seconds
- positive_rabool.
If True forces ra value to be positive
- Returns:
- rdtuple of floats
(ra, dec) - RA and Dec in degrees
- pysiaf.utils.rotations.polar_angles(vector, positive_azimuth=False)[source]¶
Compute polar coordinates of an unit vector.
- Parameters:
- vectorfloat list or array of length 3
3-component unit vector
- positive_azimuthbool
If True, the returned nu2 value is forced to be positive.
- Returns:
- nu2, nu3tuple of floats with astropy quantity
The same position represented by polar coordinates
- pysiaf.utils.rotations.posangle(attitude, v2, v3)[source]¶
Return the V3 angle at arbitrary v2,v3 using the attitude matrix.
This is the angle measured from North to V3 in an anti-clockwise direction i.e. North to East. Formulae from JWST-STScI-001550, SM-12, section 6.2. Subtract 1 from each index in the text to allow for python zero indexing.
- Parameters:
- attitude3 by 3 float array
the telescope attitude matrix
- v2float
V2 coordinate in arc-seconds
- v3float
V3 coordinate in arc-seconds
- Returns:
- pafloat
Angle in degrees - the position angle at (V2,V3)
- pysiaf.utils.rotations.radec(vector, positive_ra=False)[source]¶
Return RA and Dec in degrees corresponding to the unit vector vector.
- Parameters:
- vectora float array or list of length 3
represents a unit vector so should have unit magnitude if not, the normalization is forced within the method
- positive_rabool
indicating whether to force ra to be positive
- Returns:
- ra , dectuple of floats
RA and Dec in degrees corresponding to the unit vector vector
- pysiaf.utils.rotations.rodrigues(attitude)[source]¶
Interpret rotation matrix as a single rotation by angle phi around unit length axis.
Return axis, angle and matching quaternion. The quaternion is given in a slightly irregular order with the angle value preceding the axis information. Most of the literature shows the reverse order but JWST flight software uses the order given here.
- Parameters:
- attitude3 by 3 float array
the telescope attitude matrix
- Returns:
- axisfloat array of length 3
a unit vector which is the rotation axis
- phifloat
angle of rotation in degrees
- qfloat array of length 4
the equivalent quaternion
- pysiaf.utils.rotations.rotate(axis, angle)[source]¶
Implement fundamental 3D rotation matrices.
Rotate a vector by an angle measured in degrees, about axis 1 2 or 3 in the inertial frame. This is an anti-clockwise rotation when sighted along the axis, commonly called a right-handed rotation.
- Parameters:
- axisint
axis number, 1, 2, or 3
- anglefloat
angle of rotation in degrees
- Returns:
- rfloat array
a (3 x 3) matrix which performs the specified rotation.
- pysiaf.utils.rotations.rv(v2, v3)[source]¶
Rotate from v2,v3 position to V1 axis.
Rotate so that a V2,V3 position ends up where V1 started.
- Parameters:
- v2float
V2 position in arc-sec
- v3float
V3 position in arc-sec
- Returns:
- rva (3 x 3) array
matrix which performs the rotation described.
- pysiaf.utils.rotations.sky_posangle(attitude, ra, dec)[source]¶
Return the V3 angle at arbitrary RA and Dec using the attitude matrix.
This is the angle measured from North to V3 in an anti-clockwise direction.
- Parameters:
- attitude3 by 3 float array
the telescope attitude matrix
- rafloat
RA position in degrees
- decfloat
Dec position in degrees
- Returns:
- pafloat
resulting position angle in degrees
- pysiaf.utils.rotations.sky_to_tel(attitude, ra, dec, verbose=False)[source]¶
Transform from sky (RA, Dec) to telescope (nu2, nu3) angles.
Return nu2,nu3 position on the idealized focal sphere of any RA and Dec using the inverse of attitude matrix.
- Parameters:
- attitude3 by 3 float array
The attitude matrix.
- rafloat (default unit is degree)
RA of sky position
- decfloat (default unit is degree)
Dec of sky position
- Returns:
- nu2, nu3tuple of floats with quantity
spherical coordinates at matching position on the idealized focal sphere
- pysiaf.utils.rotations.slew(v2t, v3t, v2a, v3a)[source]¶
Calculate matrix which slews from target (v2t,v3t) to aperture position (v2a, v3a) without a roll change.
Useful for target acquisition calculations.
- Parameters:
- v2tfloat
Initial V2 position in arc-sec
- v3tfloat
Initial V3 position in arc-sec
- v2afloat
Final V2 position in arc-sec
- v3afloat
Final V3 position in arc-sec
- Returns:
- mva (3 x 3) float array
The matrix that performs the rotation described
- pysiaf.utils.rotations.tel_to_sky(attitude, nu2, nu3, positive_ra=True)[source]¶
Calculate where a nu2,nu3 position points on the sky.
- Parameters:
- attitude3 by 3 float array
the telescope attitude matrix
- nu2float or array of floats (default unit is arcsecond)
V2 coordinate in arc-seconds
- nu3float or array of floats (default unit is arcsecond)
V3 coordinate in arc-seconds
- positive_rabool.
If True forces ra value to be positive
- Returns:
- rdtuple of floats with quantity
(ra, dec) - RA and Dec
- pysiaf.utils.rotations.unit(ra, dec)[source]¶
Convert inertial frame vector expressed in polar coordinates / Euler angles to unit vector components.
See Section 5 of JWST-STScI-001550 and Equation 4.1 of JWST-PLAN-006166.
- Parameters:
- rafloat or array of floats
RA of sky position in degrees
- decfloat or array of floats
Dec of sky position in degrees
- Returns:
- vectorfloat array of length 3
the equivalent unit vector in the inertial frame
- pysiaf.utils.rotations.unit_vector_from_cartesian(x=None, y=None, z=None)[source]¶
Return unit vector corresponding to two cartesian coordinates.
Array inputs are supported.
- Parameters:
- xfloat or quantity
cartesian unit vector X coordinate in radians
- yfloat or quantity
cartesian unit vector Y coordinate in radians
- zfloat or quantity
cartesian unit vector Z coordinate in radians
- Returns:
- unit_vectornumpy.ndarray
Unit vector
- pysiaf.utils.rotations.unit_vector_hst_fgs_object(rho, phi)[source]¶
Return unit vector on the celestial sphere.
This is according to the HST object space angle definitions, CSC/TM-82/6045 1987, Section 4.1.2.2.4
- Parameters:
- rhofloat or array of floats (default unit is degree)
RA of sky position
- phifloat or array of floats (default unit is degree)
Dec of sky position
- Returns:
- vectorfloat array of length 3
the equivalent unit vector in the inertial frame
- pysiaf.utils.rotations.unit_vector_sky(ra, dec)[source]¶
Return unit vector on the celestial sphere.
- Parameters:
- rafloat or array of floats (default unit is degree)
RA of sky position
- decfloat or array of floats (default unit is degree)
Dec of sky position
- Returns:
- vectorfloat array of length 3
the equivalent unit vector in the inertial frame
- pysiaf.utils.rotations.v2v3(vector)[source]¶
Compute v2,v3 polar coordinates corresponding to a unit vector in the rotated (telescope) frame.
See Section 5 of JWST-STScI-001550.
- Parameters:
- vectorfloat list or array of length 3
unit vector of cartesian coordinates in the rotated (telescope) frame
- Returns:
- v2, v3tuple of floats
The same position represented by V2, V3 values in arc-sec.
siaf¶
Module to handle Science Instrument Aperture Files (SIAF).
The siaf module defined classes and functions to support working with SIAF files. The main class is ApertureCollection, and the Siaf class inherits from it. The class methods support basic operations and plotting.
- ApertureCollection is essentially a container for a set of pysiaf
aperture objects.
References¶
Some of the Siaf class methods were adapted from the jwxml package (https://github.com/mperrin/jwxml).
- class pysiaf.siaf.ApertureCollection(aperture_dict=None)[source]¶
Structure class for an aperture collection, e.g. read from a SIAF file.
Initialize and generate table of contents.
- class pysiaf.siaf.Siaf(instrument, filename=None, basepath=None, AperNames=None)[source]¶
Science Instrument Aperture File class.
This is a class interface to SIAF information, e.g. stored in an XML file in the PRD. It enables apertures retrieval by name, plotting, and other functionality. See the Aperture class for the detailed implementation of the transformations.
Adapted from https://github.com/mperrin/jwxml
The HST case is treated here as an instrument, because it’s single SIAF contains all apertures of all HST-instruments
Examples
fgs_siaf = SIAF(‘FGS’) fgs_siaf.apernames # returns a list of aperture names ap = fgs_siaf[‘FGS1_FULL’] # returns an aperture object ap.plot(frame=’Tel’) # plot one aperture fgs_siaf.plot() # plot all apertures in this file
- Attributes:
- observatorystr
Name of observatory
Read a SIAF from disk.
- Parameters:
- instrumentstring
one of ‘NIRCam’, ‘NIRSpec’, ‘NIRISS’, ‘MIRI’, ‘FGS’; case-insensitive.
- basepathstring
Directory to look in for SIAF files
- filenamestring, optional
Alternative method to specify a specific SIAF XML file.
- property apernames¶
List of aperture names defined in this SIAF.
- delete_aperture(aperture_name)[source]¶
Remove an aperture from the Siaf.
- Parameters:
aperture_name – str or list
- Returns:
- plot(frame='tel', names=None, label=False, units=None, clear=True, show_frame_origin=None, mark_ref=False, subarrays=True, ax=None, **kwargs)[source]¶
Plot all apertures in this SIAF.
- Parameters:
- nameslist of strings
A subset of aperture names, if you wish to plot only a subset
- subarraysbool
Plot all the minor subarrays if True, else just plot the “main” apertures
- labelbool
Add text labels stating aperture names
- unitsstr
one of ‘arcsec’, ‘arcmin’, ‘deg’
- clearbool
Clear plot before plotting (set to false to overplot)
- show_frame_originstr or list
Plot frame origin (goes to plot_frame_origin()): None, ‘all’, ‘det’, ‘sci’, ‘raw’, ‘idl’, or a list of these.
- mark_refbool
Add markers for the reference (V2Ref, V3Ref) point in each apertyre
- framestr
Which coordinate system to plot in: ‘tel’, ‘idl’, ‘sci’, ‘det’
- axmatplotlib.Axes
Desired destination axes to plot into (If None, current axes are inferred from pyplot.)
- Other matplotlib standard parameters may be passed in via **kwargs
- to adjust the style of the displayed lines.
- plot_detector_channels(frame=None, ax=None)[source]¶
Mark on the plot the various detector readout channels.
These are depicted as alternating light/dark bars to show the regions read out by each of the output amps.
- Parameters:
- framestr
Which coordinate system to plot in: ‘Tel’, ‘Idl’, ‘Sci’, ‘Det’ Optional if you have already called plot() to specify a coordinate frame.
- axmatplotlib.Axes
Desired destination axes to plot into (If None, current axes are inferred from pyplot.)
- plot_frame_origin(frame=None, which='sci', units='arcsec', ax=None)[source]¶
Mark on the plot the frame’s origin in Det and Sci coordinates.
- Parameters:
- framestr
Which coordinate system to plot in: ‘tel’, ‘idl’, ‘sci’, ‘det’ Optional if you have already called plot() to specify a coordinate frame.
- whichstr or list
Which origin to plot: ‘all’, ‘det’, ‘sci’, ‘raw’, ‘idl’, or a list
- unitsstr
one of ‘arcsec’, ‘arcmin’, ‘deg’
- axmatplotlib.Axes
Desired destination axes to plot into (If None, current axes are inferred from pyplot.)
- pysiaf.siaf.get_jwst_apertures(apertures_dict, include_oss_apertures=False, exact_pattern_match=False)[source]¶
Return ApertureCollection that corresponds to constraints specified in apertures_dict.
- Parameters:
- apertures_dictdict
Dictionary of apertures
- include_oss_aperturesbool
Whether to include OSS apertures
- exact_pattern_matchbool
- Returns:
- ApertureCollection
ApertureCollection
object Collection of apertures corresponding to selection criteria
- ApertureCollection
Examples
apertures_dict = {‘instrument’:[‘FGS’]} apertures_dict[‘pattern’] = [‘FULL’]*len(apertures_dict[‘instrument’]) fgs_apertures_all = get_jwst_apertures(apertures_dict)
- pysiaf.siaf.plot_all_apertures(subarrays=True, showorigin=True, detector_channels=True, **kwargs)[source]¶
Plot all apertures.
- pysiaf.siaf.plot_main_apertures(label=False, darkbg=False, detector_channels=False, frame='tel', attitude_matrix=None, **kwargs)[source]¶
Plot main/master apertures.
- Parameters:
- framestring
Either ‘tel’ or ‘sky’. (It does not make sense to plot apertures from multiple instruments in any of the other frames)
- attitude_matrix3x3 ndarray
Rotation matrix representing observatory attitude. Needed for sky frame plots.
tools¶
A collection of helper functions to support pysiaf calculations and functionalities.
- pysiaf.utils.tools.compute_roundtrip_error(A, B, C, D, offset_x=0.0, offset_y=0.0, verbose=False, instrument='', grid_amplitude=None)[source]¶
Return the roundtrip error of the distortion transformations specified by A,B,C,D.
Test whether the forward and inverse idl-sci transformations are consistent.
- Parameters:
- Anumpy array
polynomial coefficients
- Bnumpy array
polynomial coefficients
- Cnumpy array
polynomial coefficients
- Dnumpy array
polynomial coefficients
- offset_xfloat
Offset subtracted from input coordinate
- offset_yfloat
Offset subtracted from input coordinate
- verbosebool
verbosity
- instrumentstr
Instrument name
- Returns:
- error_estimation_metric, dx.mean(), dy.mean(), dx.std(), dy.std(), datatuple
mean and std of errors, data used in computations
- pysiaf.utils.tools.convert_polynomial_coefficients(A_in, B_in, C_in, D_in, oss=False, inverse=False, parent_aperture=None)[source]¶
Emulate some transformation made in nircam_get_polynomial_both.
Written by Johannes Sahlmann 2018-02-18, structure largely based on nircamtrans.py code by Colin Cox.
- Parameters:
- A_innumpy array
polynomial coefficients
- B_innumpy array
polynomial coefficients
- C_innumpy array
polynomial coefficients
- D_innumpy array
polynomial coefficients
- ossbool
Whether this is an OSS aperture or not
- inversebool
Whether this is forward or backward/inverse transformation
- parent_aperturestr
Name of parent aperture
- Returns:
- AR, BR, CR, DR, V3SciXAngle, V3SciYAngle, V2Ref, V3Reftuple of arrays and floats
Converted polynomial coefficients
- pysiaf.utils.tools.correct_V3SciXAngle(V3SciXAngle_deg)[source]¶
Correct input angle.
- Parameters:
- V3SciXAngle_deg
- Returns:
- V3SciXAngle_degfloat
- pysiaf.utils.tools.correct_V3SciYAngle(V3SciYAngle_deg)[source]¶
Correct input angle.
- Parameters:
- V3SciYAngle_deg
- Returns:
- V3SciYAngle_deg_correctedfloat
- pysiaf.utils.tools.get_grid_coordinates(n_side, centre, x_width, y_width=None, max_radius=None)[source]¶
Return tuple of arrays that contain the coordinates on a regular grid.
- Parameters:
- n_side: int
Number of points per side. The returned arrays have n_side**2 entries.
- centre: tuple of floats
Center coordinate
- x_width: float
Extent of the grid in the first dimension
- t_width: float
Extent of the grid in the second dimension
- Returns:
- xarray
- yarray
- pysiaf.utils.tools.is_ipython()[source]¶
Function that returns True if the user is in an ipython session and False if they are not
- pysiaf.utils.tools.jwst_fgs_to_fgs_matrix(direction='fgs2_to_fgs1', siaf=None, verbose=False)[source]¶
Return JWST FGS1_OSS to FGS2_OSS transformation matrix as stored in LoadsPRD.
- Parameters:
- siafpysiaf.Siaf instance
JWST FGS SIAF content
- Returns:
- R12ndarray
rotation matrix
References
- pysiaf.utils.tools.match_v2v3(aperture_1, aperture_2, verbose=False, match_v2_only=False)[source]¶
Modify the X[Y]DetRef,X[Y]SciRef attributes of aperture_2 such that V2Ref,V3Ref of both apertures match.
Also shift the polynomial coefficients to reflect the new reference point origin and for NIRCam recalculate angles.
- Parameters:
- aperture_1
pysiaf.Aperture object
Aperture whose V2,V3 reference position is to be used
- aperture_2
pysiaf.Aperture object
The V2,V3 reference position is to be altered to match that of aperture_1
- verbosebool
verbosity
- aperture_1
- Returns:
- new_aperture_2:
pysiaf.Aperture object
An aperture object derived from aperture_2 but with some parameters changed to match altered V2V3.
- new_aperture_2:
- pysiaf.utils.tools.revert_correct_V3SciXAngle(V3SciXAngle_deg)[source]¶
Return corrected V3SciXAngle.
- Parameters:
- V3SciXAngle_degfloat
Angle in deg
- Returns:
- V3SciXAngle_degfloat
Angle in deg
- pysiaf.utils.tools.revert_correct_V3SciYAngle(V3SciYAngle_deg)[source]¶
Return corrected V3SciYAngle.
Only correct if the original V3SciYAngle in [0,180) deg
- Parameters:
- V3SciYAngle_degfloat
angle in deg
- Returns:
- V3SciYAngle_degfloat
Angle in deg
- pysiaf.utils.tools.set_reference_point_and_distortion(instrument, aperture, parent_aperture)[source]¶
Set V2Ref and V3ref and distortion coefficients for an aperture with a parent_aperture.
- Parameters:
- instrumentstr
Instrument name
- aperture
pysiaf.Aperture
object Aperture
- parent_aperture
pysiaf.Aperture
object Parent aperture
- pysiaf.utils.tools.v3sciyangle_to_v3idlyangle(v3sciyangle)[source]¶
Convert V3SciYAngle to V3IdlYAngle.
- Parameters:
- v3sciyanglefloat
angle
- Returns:
- v3sciyanglefloat
angle
- pysiaf.utils.tools.write_matrix_to_file(matrix, file, comments=None, format='jwst_fsw_patch_request')[source]¶
Write the elements of a matrix to a text file.
- Parameters:
- matrixndarray
the matrix
- filestr
output file name
- commentsdict
comments to include in the commented file header
- formatstr
Formatting of matrix elements in output
write¶
Functions to write Science Instrument Aperture Files (SIAF).
SIAF content in an aperture_collection object can be written to an xml file that can be ingested in the PRD. Format and order of the xml fields are defined in SIAF reference files. Writing to Microsoft Excel .xlsx format is supported. Writing to .csv and other formats supported by astropy.table.Table.write is enabled.
- pysiaf.iando.write.write_jwst_siaf(aperture_collection, filename=None, basepath=None, label=None, file_format='xml', verbose=True)[source]¶
Write the content of aperture_collection into xml and xlsx files that are PRD-compliant.
- Parameters:
- aperture_collectionApertureCollection
dictionary of apertures
- filenamestr
The file name and path if you do not wish to use the default naming
- basepathstr
If you wish to use the default naming, basepath allows you to set the path where the file will be saved
- labelstr
Append default file name (“INSTR_SIAF”) with this string
- file_formatstr list
one of [‘xml’, ‘xlsx’, ‘csv’, and formats supported by astropy Table.write]
- verbose
- Returns:
- filenameslist
list of the filenames written out