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.