sorcha.ephemeris

Submodules

Attributes

AU_KM

AU_M

RADIUS_EARTH_KM

SPEED_OF_LIGHT

OBLIQUITY_ECLIPTIC

Classes

Observatory

Class containing various utility tools related to the calculation of the observatory position

Functions

create_ecl_to_eq_rotation_matrix(ecl)

Creates a rotation matrix for transforming ecliptical coordinates

make_retriever(→ pooch.Pooch)

Helper function that will create a Pooch object to track and retrieve files.

barycentricObservatoryRates(et, obsCode, observatories)

Computes the position and rate of motion for the observatory in barycentric coordinates

ecliptic_to_equatorial(v[, rot_mat])

Converts an ecliptic-aligned vector to an equatorially-aligned vector

integrate_light_time(sim, ex, t, r_obs[, lt0, iter, ...])

Performs the light travel time correction between object and observatory iteratively for the object at a given reference time

ra_dec2vec(ra, dec)

Converts a RA/Dec pair to a unit vector on the sphere

mjd_tai_to_epoch(mjd_tai)

Converts a MJD value in TAI to SPICE ephemeris time

parse_orbit_row(row, epochJD_TDB, ephem, sun_dict, ...)

Parses the input orbit row, converting it to the format expected by

create_assist_ephemeris(→ tuple)

Build the ASSIST ephemeris object

furnish_spiceypy(args, auxconfigs)

Builds the SPICE kernel, downloading the required files if needed

precompute_pointing_information(pointings_df, args, ...)

This function is meant to be run once to prime the pointings dataframe

create_ephemeris(orbits_df, pointings_df, args, sconfigs)

Generate a set of observations given a collection of orbits

universal_cartesian(mu, q, e, incl, longnode, argperi, ...)

Converts from a series of orbital elements into state vectors

universal_cometary(mu, x, y, z, vx, vy, vz, epochMJD_TDB)

Converts from a state vectors into cometary orbital elements

Package Contents

AU_KM = 149597870.7[source]
AU_M = 149597870700[source]
RADIUS_EARTH_KM = 6378.137[source]
SPEED_OF_LIGHT = 173.1446326742403[source]
OBLIQUITY_ECLIPTIC[source]
create_ecl_to_eq_rotation_matrix(ecl)[source]

Creates a rotation matrix for transforming ecliptical coordinates to equatorial coordinates. A rotation matrix based on the solar system's ecliptic obliquity is already provided as ECL_TO_EQ_ROTATION_MATRIX.

Parameters:

ecl (float) -- The ecliptical obliquity.

Returns:

rotmat -- rotation matrix for transofmring ecliptical coordinates to equatorial coordinates. Array has shape (3,3).

Return type:

numpy array/matrix of floats

make_retriever(auxconfigs, directory_path: str = None) pooch.Pooch[source]

Helper function that will create a Pooch object to track and retrieve files.

Parameters:
  • auxconfigs (dataclass) -- Dataclass of auxiliary configuration file arguments.

  • directory_path (string, default=None) -- The base directory to place all downloaded files.

Returns:

The instance of a Pooch object used to track and retrieve files.

Return type:

pooch

barycentricObservatoryRates(et, obsCode, observatories, Rearth=RADIUS_EARTH_KM, delta_et=10)[source]

Computes the position and rate of motion for the observatory in barycentric coordinates

Parameters:
  • et (float) -- JPL ephemeris time

  • obsCode (str) -- MPC observatory code

  • observatories (Observatory) -- Observatory object with spherical representations for the obsCode

  • Rearth (float, default=RADIUS_EARTH_KM) -- Radius of the Earth units[km]

  • delta_et (float, default=10) -- Difference in ephemeris time (in seconds) to derive the rotation matrix from the fixed Earth equatorial frame to J2000

Returns:

  • array -- Position of the observatory (baricentric)

  • array -- Velocity of the observatory (baricentric)

ecliptic_to_equatorial(v, rot_mat=ECL_TO_EQ_ROTATION_MATRIX)[source]

Converts an ecliptic-aligned vector to an equatorially-aligned vector

Parameters:
  • v (array (3 entries)) -- vector

  • rot_mat (2D array (3x3 matrix), default=ECL_TO_EQ_ROTATION_MATRIX) -- Rotation matrix. Default is the matrix that computes the ecliptic to equatorial conversion

Returns:

v -- Rotated vector

Return type:

array (3 entries)

integrate_light_time(sim, ex, t, r_obs, lt0=0, iter=3, speed_of_light=SPEED_OF_LIGHT)[source]

Performs the light travel time correction between object and observatory iteratively for the object at a given reference time

Parameters:
  • sim (simulation) -- Rebound simulation object

  • ex (simulation extras) -- ASSIST simulation extras

  • t (float) -- Target time

  • r_obs (array (3 entries)) -- Observatory position at time t

  • lt0 (float, default=0) -- First guess for light travel time

  • iter (int, default=3) -- Number of iterations

  • speed_of_light (float, default=SPEED_OF_LIGHT) -- Speed of light for the calculation (default is SPEED_OF_LIGHT constant)

Returns:

  • rho (array) -- Object-observatory vector

  • rho_mag (float) -- Magnitude of rho vector

  • lt (float) -- Light travel time

  • target (array) -- Object position vector at t-lt

  • vtarget (array) -- Object velocity at t-lt

ra_dec2vec(ra, dec)[source]

Converts a RA/Dec pair to a unit vector on the sphere :param ra: Target RA :type ra: float :param dec: Target dec :type dec: float

Returns:

Unit vector

Return type:

array

mjd_tai_to_epoch(mjd_tai)[source]

Converts a MJD value in TAI to SPICE ephemeris time

Parameters:

mjd_tai (float) -- Input mjd

Returns:

epoch

Return type:

Ephemeris time

class Observatory(args, auxconfigs, oc_file=None)[source]

Class containing various utility tools related to the calculation of the observatory position

observatoryPositionCache
ObservatoryXYZ
convert_to_geocentric(obs_location: dict) tuple[source]

Converts the observatory location to geocentric coordinates

Parameters:

obs_location (dict) -- Dictionary with Longitude and sin/cos of the observatory Latitude

Returns:

Geocentric position (x,y,z)

Return type:

tuple

barycentricObservatory(et, obsCode, Rearth=RADIUS_EARTH_KM)[source]

Computes the barycentric position of the observatory

Parameters:
  • et (float) -- JPL internal ephemeris time

  • obsCode (str) -- MPC Observatory code

  • Rearth (float default=RADIUS_EARTH_KM) -- Radius of the Earth, units[km]

Returns:

Barycentric position of the observatory (x,y,z)

Return type:

array (3,)

parse_orbit_row(row, epochJD_TDB, ephem, sun_dict, gm_sun, gm_total)[source]

Parses the input orbit row, converting it to the format expected by the ephemeris generation code later on

Parameters:
  • row (Pandas dataframe row) -- Row of the input dataframe

  • epochJD_TDB (float) -- epoch of the elements, in JD TDB

  • ephem (Ephem) -- ASSIST ephemeris object

  • sun_dict (dict) -- Dictionary with the position of the Sun at each epoch

  • gm_sun (float) -- Standard gravitational parameter GM for the Sun

  • gm_total (float) -- Standard gravitational parameter GM for the Solar System barycenter

Returns:

State vector (position, velocity)

Return type:

tuple

create_assist_ephemeris(args, auxconfigs) tuple[source]

Build the ASSIST ephemeris object Parameter --------- args: dictionary

Dictionary of command-line arguments.

auxconfigs: dataclass

Dataclass of auxiliary configuration file arguments.

Returns:

  • Ephem (ASSIST ephemeris object) -- The ASSIST ephemeris object

  • gm_sun (float) -- value for the GM_SUN value

  • gm_total (float) -- value for gm_total

furnish_spiceypy(args, auxconfigs)[source]

Builds the SPICE kernel, downloading the required files if needed :param args: dictionary of command-line arguments. :type args: dictionary or sorchaArguments object :param auxconfigs: Dataclass of auxiliary configuration file arguments. :type auxconfigs: dataclass

precompute_pointing_information(pointings_df, args, sconfigs)[source]

This function is meant to be run once to prime the pointings dataframe with additional information that Assist & Rebound needs for it's work.

Parameters:
  • pointings_df (pandas dataframe) -- Contains the telescope pointing database.

  • args (dictionary) -- Command line arguments needed for initialization.

  • sconfigs (dataclass) -- Dataclass of configuration file arguments.

Returns:

pointings_df -- The original dataframe with several additional columns of precomputed values.

Return type:

pandas dataframe

create_ephemeris(orbits_df, pointings_df, args, sconfigs)[source]

Generate a set of observations given a collection of orbits and set of pointings.

Parameters:
  • orbits_df (pandas dataframe) -- The dataframe containing the collection of orbits.

  • pointings_df (pandas dataframe) -- The dataframe containing the collection of telescope/camera pointings.

  • args -- Various arguments necessary for the calculation

  • sconfigs --

    Dataclass of configuration file arguments. Various configuration parameters necessary for the calculation ang_fov : float

    The angular size (deg) of the field of view

    bufferfloat

    The angular size (deg) of the buffer around the field of view. A buffer is required to allow for some motion between the time of the observation and the time of the picket (t_picket)

    picket_intervalfloat

    The interval (days) between picket calculations. This is 1 day by default. Current there is only one such interval, used for all objects. It is currently possible for extremely fast-moving objects to be missed. This will be remedied in future releases.

    obsCodestring

    The MPC code for the observatory. (This is current a configuration parameter, but these should be included in the visit information, to allow for multiple observatories.

    nsideinteger

    The nside value used for the HEALPIx calculations. Must be a power of 2 (1, 2, 4, ...) nside=64 is current default.

    n_sub_intervals: int

    Number of sub-intervals for the Lagrange interpolation (default: 101)

Returns:

observations -- The dataframe of observations needed for Sorcha to continue

Return type:

pandas dataframe

Notes

This works by calculating and regularly updating the sky-plane locations (unit vectors) of all the objects in the collection of orbits. The HEALPix index for each of the locations is calculated. A dictionary with pixel indices as keys and lists of ObjIDs for those objects in each HEALPix tile as values is generated. An individual one of these calculations is called a 'picket', as one element of a long picket fence. Typically, the interval between pickets is one day.

Given a specific pointing, the set of HEALPix tiles that are overlapped by the pointing (and a buffer region) is computed. Then the precise locations of just those objects within that set of HEALPix tiles are computed. Details for those that actually do land within the field of view are passed along.

universal_cartesian(mu, q, e, incl, longnode, argperi, tp, epochMJD_TDB)[source]

Converts from a series of orbital elements into state vectors using the universal variable formulation

The output vector will be oriented in the same system as the positional angles (i, Omega, omega)

Note that mu, q, tp and epochMJD_TDB must have compatible units As an example, if q is in au and tp/epoch are in days, mu must be in (au^3)/days^2

Parameters:
  • mu (float) -- Standard gravitational parameter GM (see note above about units)

  • q (float) -- Perihelion (see note above about units)

  • e (float) -- Eccentricity

  • incl (float) -- Inclination (radians)

  • longnode (float) -- Longitude of ascending node (radians)

  • argperi (float) -- Argument of perihelion (radians)

  • tp (float) -- Time of perihelion passage in TDB scale (see note above about units)

  • epochMJD_TDB (float) -- Epoch (in TDB) when the elements are defined (see note above about units)

Returns:

  • float -- x coordinate

  • float -- y coordinate

  • float -- z coordinate

  • float -- x velocity

  • float -- y velocity

  • float -- z velocity

universal_cometary(mu, x, y, z, vx, vy, vz, epochMJD_TDB)[source]

Converts from a state vectors into cometary orbital elements using the universal variable formulation

The input vector will determine the orientation of the positional angles (i, Omega, omega)

Note that mu and the state vectors must have compatible units As an example, if x is in au and vx are in au/days, mu must be in (au^3)/days^2

Parameters:
  • mu (float) -- Standard gravitational parameter GM (see note above about units)

  • x (float) -- x coordinate

  • y (float) -- y coordinate

  • z (float) -- z coordinate

  • vx (float) -- x velocity

  • vy (float) -- y velocity

  • vz (float) -- z velocity

  • (float) (epochMJD_TDB) -- Epoch (in TDB) when the elements are defined (see note above about units)

Returns:

  • q (float) -- Perihelion (see note above about units)

  • e (float) -- Eccentricity

  • incl (float) -- Inclination (radians)

  • longnode (float) -- Longitude of ascending node (radians)

  • argperi (float) -- Argument of perihelion (radians)

  • tp (float) -- Time of perihelion passage in TDB scale (see note above about units)