rss_ringoccs.occgeo.calc_occ_geometry module¶
Purpose: | Functions for calculating occultation geometry parameters listed in *GEO.LBL file from Archived_Cassini_RSS_RingOccs_2018/ and other useful geometry information, such as free-space regions and planetary occultation times. |
---|---|
Dependencies: |
|
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_B_deg
(et_vals, spacecraft, dsn, nhat_p, kernels=None, ref='J2000')¶ This calculates ring opening angle, or the observed ring elevation, as the complement to the angle made by the planet pole vector and the spacecraft to DSN vector
- Arguments
et_vals (np.ndarray): Array of earth-received times in ephemeris seconds. dsn (str): DSN observing station ID – must be compatible with NAIF. nhat_p (np.ndarray): 1x3 array unit vector in planet pole direction. - Keyword Arguments
kernels (str or list): List of NAIF kernels, including path. ref (str): Reference frame to be used in spiceypy calls. Default is ‘J2000’ - Returns
B_deg_vals (np.ndarray): Array of ring opening angle in degrees.
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_B_eff_deg
(B_deg, phi_ora_deg)¶ This calculates the effective ring opening angle in degrees.
- Arguments
B_deg (float or np.ndarray): Ring opening angle in degrees. phi_ora_deg (float or np.ndarray): Observed ring azimuth in degrees. - Returns
B_eff_deg (float or np.ndarray): Effective ring opening angle in degrees. - Notes
- Reference: [GRESH86] Eq. 16
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_D_km
(t1, t2)¶ This calculates the light distance between two input times, in km.
- Args
t1 (np.ndarray): Array of time in seconds. t2 (np.ndarray): Array of time in seconds. - Returns
D_km_vals (np.ndarray): Array of light distance in km.
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_F_km
(D_km_vals, f_sky_hz_vals, B_deg_vals, phi_ora_deg_vals)¶ This calculates the Fresnel scale using Eq. 6 of [MTR1986].
- Arguments
D_km_vals (np.ndarray): Array of spacecraft to ring intercept point distances in km. f_sky_hz_vals (np.ndarray): Array of downlink sinusoidal signal frequency at the front-end of observing dsn station, in Hz. B_deg_vals (np.ndarray): Array of ring opening angle in degrees. phi_ora_deg_vals (np.ndarray): Array of observed ring azimuth in degrees. - Returns
F_km_vals (np.ndarray): Array of Fresnel scale in km. - Notes
- diffcorr uses an independently-calculated Fresnel scale
- Reference: [MTR1986] Equation 6
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_beta
(B_deg, phi_ora_deg)¶ This calculates the optical depth enhancement factor.
- Arguments
B_deg (float or np.ndarray): Ring opening angle in degrees. phi_ora_deg (float or np.ndarray): Observed ring azimuth in degrees. - Returns
beta (np.ndarray): Optical depth enhancement factor.
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_elevation_deg
(et_vals, target, obs, kernels=None)¶ Calculate the elevation of a target above the horizon for a given observer.
- Arguments
et_vals (np.ndarray): Array of observed event times in ET sec. target (str): Target name – must be compatible with NAIF. This will typically be spacecraft or planet name. obs (str): Observer name – must be compatible with NAIF. This will typically be observing dsn station name. Observer is assumed to be on Earth. kernels (str or list): List of NAIF kernels, including path. - Returns
elev_deg_vals (np.ndarray): Array of elevation angles in degrees.
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_impact_radius_km
(R_sc_km_vals, et_vals, spacecraft, dsn, nhat_p, ref='J2000', kernels=None)¶ This calculates the closest approach of the spacecraft signal to the planet defined as a sphere.
- Arguments
R_sc_km_vals (list): List of 3-element arrays of spacecraft position vector in planetocentric frame at input et_vals. et_vals (np.ndarray): Array of Earth-received times in ephemeris seconds. spacecraft (str): Spacecraft name dsn (str): DSN observing station ID nhat_p (np.ndarray): 1x3 array unit vector in planet pole direction. - Keyword Arguments
kernels (list): List of NAIF kernels, including path. ref (str): Reference frame to be used in spiceypy calls. Default is ‘J2000’ - Returns
R_imp_km_vals (np.ndarray): Array of impact radius in km.
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_phi_deg
(et_vals, rho_vec_km_vals, spacecraft, dsn, nhat_p, ref='J2000', kernels=None)¶ This calculates observed ring azimuth and ring longitude.
- Arguments
et_vals (np.ndarray): Array of earth-received time in ET seconds. rho_vec_km_vals (np.ndarray): Nx3 array of ring intercept position vectors in km. spacecraft (str): Name of spacecraft dsn (str): DSN observing station ID nhat_p (np.ndarray): 1x3 array unit vector in planet pole direction. - Keyword Arguments
kernels (str or list): List of NAIF kernels, including path ref (str): Reference frame to be used in spiceypy calls. Default is ‘J2000’ - Returns
phi_rl_deg_vals (np.ndarray): Array of inertial longitude in degrees. phi_ora_deg_vals (np.ndarray): Array of observed ring azimuth in degrees. - Notes:
- phi_ora_deg differs from the [MTR1986] definition by 180 degrees.
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_rho_km
(et_vals, planet, spacecraft, dsn, kernels=None, ring_frame=None)¶ Calculate the distance between Saturn center to ring intercept point.
- Arguments
et_vals (np.ndarray): Array of observed event times in ET sec. planet (str): Planet name spacecraft (str): Spacecraft name dsn (str): DSN observing station ID - Keyword Arguments
kernels (str or list): List of NAIF kernels, including path. ring_frame (str): Ring plane frame. Default is the equatorial frame, (e.g. ‘IAU_SATURN’) - Returns
rho_km_vals (np.ndarray): Array of ring intercept points in km.
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_rho_vec_km
(et_vals, planet, spacecraft, dsn, ref='J2000', kernels=None, verbose=False, ring_frame=None)¶ This calculates the position vector of the ring intercept point from the planet center in J2000 frame.
- Arguments
et_vals (np.ndarray): Array of earth-received times in ET sec planet (str): Name of planet spacecraft (str): Name of spacecraft dsn (str): DSN observing station ID - Keyword Arguments
kernels (str or list): Path to NAIF kernels verbose (bool): Option for printing processing steps ring_frame (str): Ring plane frame. Default is the equatorial frame, (e.g. ‘IAU_SATURN’) ref (str): Reference frame to be used in spiceypy calls. Default is ‘J2000’ - Returns
rho_vec_km_vals (list): List of 3xN np.ndarrays of the planet center to ring intercept point position vector in J2000 frame t_ret_et_vals (np.ndarray): Array of ring event times in ET seconds. - References
- Ring intercept point calculation using a dynamical frame. See [NAIF] page 19.
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_rip_velocity
(rho_km_vals, phi_rl_deg_vals, dt)¶ This calculates the ring intercept point radial and azimuthal velocity.
- Arguments
rho_km_vals (np.ndarray): Array of ring intercept points in km. phi_rl_deg_vals (np.ndarray): Array of ring longitudes in degrees. dt (float): Constant time spacing between points. - Returns
rho_dot_kms_vals (np.ndarray): Array of ring intercept radial velocities in km/s. phi_rl_dot_kms_vals (np.ndarray): Array of ring intercept azimuthal velocties in km/s.
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_sc_state
(et_vals, spacecraft, planet, dsn, nhat_p, ref='J2000', kernels=None)¶ This calculates spacecraft state vector in a planetocentric frame.
- Arguments
et_vals (np.ndarray): Array of spacecraft event times in ET seconds. spacecraft (str): Spacecraft name planet (str): Planet name dsn (str): Deep Space Network observing station ID nhat_p (np.ndarray): 1x3 array unit vector in planet pole direction. - Keyword Arguments
kernels (str or list): Path to NAIF kernel(s) ref (str): Reference frame to be used in spiceypy calls. Default is ‘J2000’ - Returns
R_sc_km_vals (list): List of Nx3 np.ndarrays of spacecraft position vector in km in planetocentric frame R_sc_dot_kms_vals (list): List of Nx3 np.ndarrays of spacecraft velocity vector in km/s. - Notes
- Saturn planetocentric frame is defined by x-axis parallel to projection of spacecraft-to-Earth line-of-sight, z-axis in direction of Saturn’s pole.
-
rss_ringoccs.occgeo.calc_occ_geometry.
calc_set_et
(et_vals, spacecraft, dsn, kernels=None)¶ This calculates the time at which photons left the spacecraft, in ET sec.
- Arguments
et_vals (np.ndarray): Array of earth-received times in ET seconds. spacecraft (str): dsn (str): Deep Space Network observing station ID - Keyword Arguments
kernels (str or *list): Path to NAIF kernels - Returns
t_set_et_vals (np.ndarray): Array of spacecraft event times in ET sec.
-
rss_ringoccs.occgeo.calc_occ_geometry.
find_gaps
(t_ret_spm_vals, year, doy, rho_km_vals, phi_rl_deg_vals, niter=100, tolerance=0.001, t0=2454467.0, gaps_file='../tables/gap_orbital_elements.txt', kernels=None)¶ Find regions of free-space power (gaps) in the ring system.
- Arguments
t_ret_spm_vals (np.ndarray): Ring event times in SPM year (str): Reference year for seconds past midnight doy (str): Reference day of year for seconds past midnight rho_km_vals (np.ndarray): Ring intercept points in km phi_rl_deg_vals (np.ndarray): Inertial ring longitude in deg. - Keyword Arguments
niter (int): Maximum number of iterations tolerance (float): Minimum difference between new and old guess for converging on a solution t0 (float): Reference epoch UTC 2008 January 1 12:00:00 for values in gaps_file, in Julian date gaps_file (str): Path to text file with orbital elements of Saturn ring features kernels (str or list): Path to NAIF kernels - Returns
gap_bounds (list): List of 1x2 lists of gap boundaries in km - Notes
- Reference: [NICH14]
- Given the default “gaps_file” keyword argument, this script must be run in a directory one level below the top-level rss_ringoccs directory.
-
rss_ringoccs.occgeo.calc_occ_geometry.
get_freespace
(t_ret_spm_vals, year, doy, rho_km_vals, phi_rl_deg_vals, t_oet_spm_vals, atmos_occ_spm_vals, kernels=None, split_ind=None)¶ Return list of gap boundaries (inner and outer edge) in distance from center of Saturn and in seconds past midnight.
- Arguments
t_ret_spm_vals (np.ndarray): Ring event times in SPM year (str): Reference year for seconds past midnight doy (str): Reference day of year for seconds past midnight rho_km_vals (np.ndarray): Ring intercept points in km phi_rl_deg_vals (np.ndarray): Inertial ring longitude in deg. t_oet_spm_vals (np.ndarray): Observed event times in SPM atmos_occ_spm_vals (np.ndarray): SPM times of when spacecraft signal is blocked by planet atmosphere - Keyword Arguments
kernels (str or list): Path to NAIF kernels split_ind (int): Index of when a chord event switches from ingress to egress - Returns
gaps_km (list): List of 1x2 lists of gap boundaries in km gaps_spm (list): List of 1x2 lists of gap boundaries in SPM
-
rss_ringoccs.occgeo.calc_occ_geometry.
get_freespace_km
(ret_spm, year, doy, rho_km, phi_rl_deg)¶ Get all free-space regions, in and outside ring system.
- Arguments
ret_spm (np.ndarray): Ring event times in SPM year (str): Reference year :doy (str) Reference day of year :rho_km (np.ndarray): Ring intercept points in km :phi_rl_deg (np.ndarray): Inertial ring longitudes in deg
- Returns
freespace_km (list): List of free-space boundaries in km
-
rss_ringoccs.occgeo.calc_occ_geometry.
get_planet_occ_times
(et_vals, obs, planet, spacecraft, height_above=500.0, kernels=None)¶ Return times when the spacecraft-to-observer ray is blocked by planet.
- Arguments
et_vals (np.ndarray): Array of observed event times in ET sec. obs (str): Observer name planet (str): Planet name spacecraft (str): Spacecraft name - Keyword Arguments
height_above (float): Height in km to be added to planet radius to account for the atmosphere kernels (str or list): Path to NAIF kernels - Returns
et_blocked_vals (np.ndarray): Array of observed event times in ET - Note:
- This was made to be generalizable to different planets, but has been only tested with planet=’Saturn’.
-
rss_ringoccs.occgeo.calc_occ_geometry.
get_pole
(et, planet, kernels=None)¶ Calculate unit vector in pole direction from kernel constants.
- Arguments
et (float): Ephemeris seconds past J2000 planet (str): Planet name - Keyword Arguments
kernels (str or list): Path to NAIF kernels - Returns
nhat_p (np.ndarray): 1x3 unit vector in pole direction. - Note:
- Quadratic terms for pole direction are typically zero but are retained here for consistency with PCK file format definitions.
-
rss_ringoccs.occgeo.calc_occ_geometry.
get_start_jd
(year, doy)¶ Get the start of a day in Julian date times.
- Arguments
year (str): Year doy (str): Day of year - Returns
start_jd (float): Julian date time of the start of a day
-
rss_ringoccs.occgeo.calc_occ_geometry.
rad_converge
(t_ret_spm_vals, rho_km_vals, phi_rl_deg_vals, semimajor, eccentricity, curlypi_0, curlypi_dot, niter=100, tolerance=0.001)¶ Computes initial guess for radius of ring feature using the semimajor axis. Selects time and longitude closest to guess and computes true anomaly for a new radius guess. Continues this estimation method iteratively until difference between new and old radius guesses is less than some tolerance or maximum number of iterations is reached.
- Arguments
t_ret_spm_vals (np.ndarray): Ring event times in SPM. rho_km_vals (np.ndarray): Ring intercept points in km. phi_rl_deg_vals (np.ndarray): Inertial ring longitude in deg. semimajor (float): Semimajor axis of ring feature in km. eccentricity (float): Eccentricity of ring feature. curlypi_0 (*float): Longitude of periapse in degrees. curlypi_dot (*float): Apsidal precession rate in degrees/day.
- Keyword Arguments
niter (int): Maximum number of iterations tolerance (float): Minimum difference between new and old guess for converging on a solution - Returns
radius_new (float): Estimated radius of ring feature in km - Notes
- Reference: [NICH14]
-
rss_ringoccs.occgeo.calc_occ_geometry.
remove_blocked
(t_oet_spm_vals, atmos_occ_spm_vals, t_ret_spm_vals, phi_rl_deg_vals, rho_km_vals)¶ Remove values that occur during times blocked by planet atmosphere.
- Arguments
t_oet_spm_vals (np.ndarray): Observed event times in SPM atmos_occ_spm_vals (np.ndarray): SPM times of when spacecraft signal is blocked by planet atmosphere t_ret_spm_vals (np.ndarray): Ring event times in SPM phi_rl_deg_vals (np.ndarray): Inertial ring longitude in deg. rho_km_vals (np.ndarray): Ring intercept points in km - Returns
t_ret_spm_vals (np.ndarray): Ring event times in SPM, excluding atmospheric occultation times t_oet_spm_vals (np.ndarray): Observed event times in SPM, excluding atmospheric occultation times phi_rl_deg_vals (np.ndarray): Inertial ring longitude in deg, excluding atmospheric occultation times rho_km_vals (np.ndarray): Ring intercept points in km, excluding atmospheric occultation times
-
rss_ringoccs.occgeo.calc_occ_geometry.
split_chord_arr
(t_ret_spm_vals, t_oet_spm_vals, atmos_occ_spm_vals, phi_rl_deg_vals, rho_km_vals, ind, profdir)¶ Return array of only ingress or egress portion of a chord occultation.
- Arguments
t_ret_spm_vals (np.ndarray): Ring event times in SPM t_oet_spm_vals (np.ndarray): Observed event times in SPM atmos_occ_spm_vals (np.ndarray): SPM times of when spacecraft signal is blocked by planet atmosphere phi_rl_deg_vals (np.ndarray): Inertial ring longitude in deg. rho_km_vals (np.ndarray): Ring intercept points in km ind (int): Index of where ingress switches to egress profdir (str): Profile direction to return (‘“INGRESS”’ or ‘“EGRESS”’) - Returns
t_ret_spm_vals (np.ndarray): Ring event times in SPM of ‘profdir’ portion of occultation t_oet_spm_vals (np.ndarray): Observed event times in SPM of ‘profdir’ portion of occultation phi_rl_deg_vals (np.ndarray): Inertial ring longitude in deg of ‘profdir’ portion of occultation rho_km_vals (np.ndarray): Ring intercept points in km of ‘profdir’ portion of occultation
-
rss_ringoccs.occgeo.calc_occ_geometry.
xform_j2k_to_pcf
(vec, et, spacecraft, dsn, nhat_p, ref='J2000', kernels=None)¶ Transform vector in J2000 frame to planet ring plane frame.
- Arguments
vec (np.ndarray): 3-element vector in J2000 frame et (float): ET in seconds corresponding to input vec dsn (str): DSN observing station ID nhat_p (np.ndarray): 1x3 array unit vector in planet pole direction. - Keyword Arguments
kernels (str or list): Path to NAIF kernels - Returns
out_vec (np.ndarray): 3-element vector in planet ring plane frame.