Moving Objects API

class rubin_sim.moving_objects.BaseObs(footprint='camera', r_fov=1.75, x_tol=5, y_tol=3, eph_mode='nbody', eph_type='basic', obs_code='I11', eph_file=None, obs_time_col='observationStartMJD', obs_time_scale='TAI', seeing_col='seeingFwhmGeom', visit_exp_time_col='visitExposureTime', obs_ra='fieldRA', obs_dec='fieldDec', obs_rot_sky_pos='rotSkyPos', obs_degrees=True, outfile_name='lsst_obs.dat', obs_info='', camera_footprint_file=None)[source]

Bases: object

Base class to generate observations of a set of moving objects.

Parameters:
  • footPrint (str, optional) – Specify the footprint for the FOV. Options include “camera”, “circle”, “rectangle”. ‘Camera’ means use the actual LSST camera footprint (following a rough cut with a circular FOV). Default is camera FOV.

  • r_fov (float, optional) – If footprint is “circular”, this is the radius of the fov (in degrees). Default 1.75 degrees (only used for circular fov).

  • x_tol (float, optional) – If footprint is “rectangular”, this is half of the width of the (on-sky) fov in the RA direction (in degrees). Default 5 degrees.

  • y_tol (float, optional) – If footprint is “rectangular”, this is half of the width of the fov in Declination (in degrees). Default is 3 degrees

  • eph_mode (str, optional) – Mode for ephemeris generation - nbody or 2body. Default is nbody.

  • eph_type (str, optional) – Type of ephemerides to generate - full or basic. Full includes all values calculated by openorb; Basic includes a more basic set. Default is Basic.

  • eph_file (str or None, optional) – The name of the planetary ephemerides file to use for ephemeris generation. Default (None) will use the default for PyOrbEphemerides.

  • obs_code (str, optional) – Observatory code for ephemeris generation. Default is “I11” - Cerro Pachon.

  • obs_time_col (str, optional) – Name of the time column in the obsData. Default ‘observationStartMJD’.

  • obs_time_scale (str, optional) – Type of timescale for MJD (TAI or UTC currently). Default TAI.

  • seeing_col (str, optional) – Name of the seeing column in the obsData. Default ‘seeingFwhmGeom’. This should be the geometric/physical seeing as it is used for the trailing loss calculation.

  • visit_exp_time_col (str, optional) – Name of the visit exposure time column in the obsData. Default ‘visitExposureTime’.

  • obs_ra (str, optional) – Name of the RA column in the obsData. Default ‘fieldRA’.

  • obs_dec (str, optional) – Name of the Dec column in the obsData. Default ‘fieldDec’.

  • obs_rot_sky_pos (str, optional) – Name of the Rotator column in the obsData. Default ‘rotSkyPos’.

  • obs_degrees (bool, optional) – Whether the observational data is in degrees or radians. Default True (degrees).

  • outfile_name (str, optional) – The output file name. Default is ‘lsst_obs.dat’.

  • obs_info (str, optional) – A string that captures provenance information about the observations. For example: ‘baseline_v2.0_10yrs, years 0-5’ or ‘baseline2018a minus NES’ Default ‘’.

calc_colors(sedname='C.dat', sed_dir=None)[source]

Calculate the colors for a given SED.

If the sedname is not already in the dictionary self.colors, this reads the SED from disk and calculates all V-[filter] colors for all filters in self.filterlist. The result is stored in self.colors[sedname][filter], so will not be recalculated if the SED + color is reused for another object.

Parameters:
  • sedname (str, optional) – Name of the SED. Default ‘C.dat’.

  • sed_dir (str, optional) – Directory containing the SEDs of the moving objects. Default None = $RUBIN_SIM_DATA_DIR/movingObjects,

Returns:

colors – Dictionary of the colors in self.filterlist.

Return type:

dict {‘filter’: color}}

calc_trailing_losses(velocity, seeing, texp=30.0)[source]

Calculate the detection and SNR trailing losses.

‘Trailing’ losses = loss in sensitivity due to the photons from the source being spread over more pixels; thus more sky background is included when calculating the flux from the object and thus the SNR is lower than for an equivalent brightness stationary/PSF-like source. dmagTrail represents this loss.

‘Detection’ trailing losses = loss in sensitivity due to the photons from the source being spread over more pixels, in a non-stellar-PSF way, while source detection is (typically) done using a stellar PSF filter and 5-sigma cutoff values based on assuming peaks from stellar PSF’s above the background; thus the SNR is lower than for an equivalent brightness stationary/PSF-like source (and by a greater factor than just the simple SNR trailing loss above). dmag_detect represents this loss.

Parameters:
  • velocity (np.ndarray or float) – The velocity of the moving objects, in deg/day.

  • seeing (np.ndarray or float) – The seeing of the images, in arcseconds.

  • texp (np.ndarray or float, optional) – The exposure time of the images, in seconds. Default 30.

Returns:

  • dmag Trail, dmag_detect ((np.ndarray np.ndarray))

  • or (float, float) – dmag_trail and dmag_detect for each set of velocity/seeing/texp values.

generate_ephemerides(sso, times, eph_mode=None, eph_type=None)[source]

Generate ephemerides for ‘sso’ at times ‘times’ (assuming MJDs, with timescale self.obs_time_scale).

The default engine here is pyoorb, however other ephemeris generation could be used with a matching API to PyOrbEphemerides.

The initialized pyoorb class (PyOrbEphemerides) is saved, to skip setup on subsequent calls.

Parameters:
  • sso (rubin_sim.movingObjects.Orbits) – Typically this will be a single object.

  • times (np.ndarray) – The times at which to generate ephemerides. MJD.

  • eph_mode (str or None, optional) – Potentially override default eph_mode (self.eph_mode). Must be ‘2body’ or ‘nbody’.

Returns:

ephs – Results from propigating the orbit(s) to the specified times. Columns like: obj_id, sedname, time, ra, dec, dradt, ddecdt, phase, solarelon.

Return type:

pd.Dataframe

read_filters(filter_dir=None, bandpass_root='total_', bandpass_suffix='.dat', filterlist=('u', 'g', 'r', 'i', 'z', 'y'), v_dir=None, v_filter='harris_V.dat')[source]

Read (LSST) and Harris (V) filter throughput curves.

Only the defaults are LSST specific; this can easily be adapted for any survey.

Parameters:
  • filter_dir (str, optional) – Directory containing the filter throughput curves (‘total*.dat’) Default set by ‘LSST_THROUGHPUTS_BASELINE’ env variable.

  • bandpass_root (str, optional) – Rootname of the throughput curves in filterlist. E.g. throughput curve names are bandpass_root + filterlist[i] + bandpass_suffix Default total_ (appropriate for LSST throughput repo).

  • bandpass_suffix (str, optional) – Suffix for the throughput curves in filterlist. Default ‘.dat’ (appropriate for LSST throughput repo).

  • filterlist (list, optional) – List containing the filter names to use to calculate colors. Default (‘u’, ‘g’, ‘r’, ‘i’, ‘z’, ‘y’)

  • v_dir (str, optional) – Directory containing the V band throughput curve. Default None = $RUBIN_SIM_DATA_DIR/movingObjects

  • v_filter (str, optional) – Name of the V band filter curve. Default harris_V.dat.

setup_ephemerides()[source]

Initialize the ephemeris generator. Save the setup PyOrbEphemeris class.

This uses the default engine, pyoorb - however this could be overwritten to use another generator.

sso_in_camera_fov(ephems, obs_data)[source]

Determine which observations are within the actual camera footprint for a series of observations. Note that ephems and obs_data must be the same length.

Parameters:
  • ephems (np.ndarray) – Ephemerides for the objects.

  • obs_data (np.ndarray) – Observation pointings.

Returns:

indices – Returns the indexes of the numpy array of the object observations which are inside the fov.

Return type:

np.ndarray

sso_in_circle_fov(ephems, obs_data)[source]

Determine which observations are within a circular fov for a series of observations. Note that ephems and obs_data must be the same length.

Parameters:
  • ephems (np.recarray) – Ephemerides for the objects.

  • obs_data (np.recarray) – The observation pointings.

Returns:

indices – Returns the indexes of the numpy array of the object observations which are inside the fov.

Return type:

np.ndarray

sso_in_fov(ephems, obs_data)[source]

Convenience layer - determine which footprint method to apply (from self.footprint) and use it.

Parameters:
  • ephems (np.ndarray) – Ephemerides for the objects.

  • obs_data (np.ndarray) – Observation pointings.

Returns:

indices – Returns the indexes of the numpy array of the object observations which are inside the fov.

Return type:

np.ndarray

sso_in_rectangle_fov(ephems, obs_data)[source]

Determine which observations are within a rectangular FoV for a series of observations. Note that ephems and obs_data must be the same length.

Parameters:
  • ephems (np.recarray) – Ephemerides for the objects.

  • obs_data (np.recarray) – The observation pointings.

Returns:

indices – Returns the indexes of the numpy array of the object observations which are inside the fov.

Return type:

np.ndarray

class rubin_sim.moving_objects.ChebyFits(orbits_obj, t_start, t_span, time_scale='TAI', obscode='I11', sky_tolerance=2.5, n_coeff_position=14, n_coeff_vmag=9, n_coeff_delta=5, n_coeff_elongation=6, ngran=64, eph_file=None, n_decimal=10)[source]

Bases: object

Generates chebyshev coefficients for a provided set of orbits.

Calculates true ephemerides using PyEphemerides, then fits these positions with a constrained Chebyshev Polynomial, using the routines in chebyshevUtils.py.

Parameters:
  • orbits_obj (rubin_sim.moving_objects.Orbits) – The orbits for which to fit chebyshev polynomial coefficients.

  • t_start (float) – The starting point in time to fit coefficients. MJD.

  • t_span (float) – The time span (starting at t_start) over which to fit coefficients (Days).

  • time_scale (str, optional) – One of {‘TAI’, ‘UTC’, ‘TT’} The timescale of the MJD time, t_start, and the time_scale that should be used with the chebyshev coefficients.

  • obsCode (int, optional) – The observatory code of the location for which to generate ephemerides. Default I11 (Cerro Pachon).

  • sky_tolerance (float, optional) – The desired tolerance in mas between ephemerides calculated by OpenOrb and fitted values. Default 2.5 mas.

  • nCoeff_position (int, optional) – The number of Chebyshev coefficients to fit for the RA/Dec positions. Default 14.

  • nCoeff_vmag (int, optional) – The number of Chebyshev coefficients to fit for the V magnitude values. Default 9.

  • nCoeff_delta (int, optional) – The number of Chebyshev coefficients to fit for the distance between Earth/Object. Default 5.

  • nCoeff_elongation (int, optional) – The number of Chebyshev coefficients to fit for the solar elongation. Default 5.

  • ngran (int, optional) – The number of ephemeris points within each Chebyshev polynomial segment. Default 64.

  • eph_file (str, optional) – The path to the JPL ephemeris file to use. Default is ‘$OORB_DATA/de405.dat’.

  • n_decimal (int, optional) – The number of decimal places to allow in the segment length (and thus the times of the endpoints) can be limited to n_decimal places. Default 10. For LSST SIMS moving object database, this should be 13 decimal places for NEOs and 0 for all others.

Notes

Many chebyshev polynomials are used to fit one moving object over a given timeperiod; typically, the length of each segment is typically about 2 days for MBAs. The start and end of each segment must match exactly, and the entire segments must fit into the total timespan an integer number of times. This is accomplished by setting n_decimal to the number of decimal places desired in the ‘time’ value. For faster moving objects, this number needs be greater to allow for smaller subdivisions. It’s tempting to allow flexibility to the point of not enforcing this non-overlap; however, then the resulting ephemeris may have multiple values depending on which polynomial segment was used to calculate the ephemeris.

The length of each chebyshev polynomial is related to the number of ephemeris positions used to fit that polynomial by ngran: length = timestep * ngran The length of each polynomial is adjusted so that the residuals in RA/Dec position are less than sky_tolerance - default = 2.5mas. The polynomial length (and the resulting residuals) is affected by ngran (i.e. timestep).

Default values are based on Yusra AlSayaad’s work.

calc_one_segment(orbit_obj, ephs)[source]

Calculate the coefficients for a single Chebyshev segment, for a single object.

Calculates the coefficients and residuals, and saves this information to self.coeffs, self.resids, and (if there are problems), self.failed.

Parameters:
  • orbit_obj (rubin_sim.moving_objects.Orbits) – The single Orbits object we’re fitting at the moment.

  • ephs (np.ndarray) – The ephemerides we’re fitting at the moment (for the single object / single segment).

calc_segment_length(length=None)[source]

Set the typical initial ephemeris timestep and segment length for all objects between t_start/t_end.

Sets self.length.

The segment length will fit into the time period between t_start/t_end an approximately integer multiple of times, and will only have a given number of decimal places.

Parameters:

length (float, optional) – If specified, this value for the length is used, instead of calculating it here.

calc_segments()[source]

Run the calculation of all segments over the entire time span.

generate_ephemerides(times, by_object=True)[source]

Generate ephemerides using OpenOrb for all orbits.

Parameters:

times (np.ndarray) – The times to use for ephemeris generation.

make_all_times()[source]

Using t_start and t_end, generate a numpy array containing times spaced at timestep = self.length/self.ngran. The expected use for this time array would be to generate ephemerides at each timestep.

Returns:

times – Numpy array of times.

Return type:

np.ndarray

write(coeff_file, resid_file, failed_file, append=False)[source]

Write coefficients, residuals and failed fits to disk.

Parameters:
  • coeff_file (str) – The filename for the coefficient values.

  • resid_file (str) – The filename for the residual values.

  • failed_file (str) – The filename to write the failed fit information (if failed objects exist).

  • append (bool, optional) – Flag to append (or overwrite) the output files.

class rubin_sim.moving_objects.ChebyValues[source]

Bases: object

Calculates positions, velocities, deltas, vmags and elongations, given a series of coefficients generated by ChebyFits.

get_ephemerides(times, obj_ids=None, extrapolate=False)[source]

Find the ephemeris information for ‘obj_ids’ at ‘time’.

Implicit in how this is currently written is that the segments are all expected to cover the same start/end time range across all objects. They do not have to have the same segment length for all objects.

Parameters:
  • times (float or np.ndarray) – The time to calculate ephemeris positions.

  • obj_ids (np.ndarray, opt) – The object ids for which to generate ephemerides. If None, then just uses all objects.

  • extrapolate (bool, opt) – If True, extrapolate beyond ends of segments if time outside of segment range. If False, return ValueError if time is beyond range of segments.

Returns:

ephemerides – The ephemeris positions for all objects. Note that these may not be sorted in the same order as obj_ids.

Return type:

np.ndarray

read_coefficients(cheby_fits_file)[source]

Read coefficients from output file written by ChebyFits.

Parameters:

cheby_fits_file (str) – The filename of the coefficients file.

set_coefficients(cheby_fits)[source]

Set coefficients using a ChebyFits object. (which contains a dictionary of obj_id, t_start, t_end, ra, dec, delta, vmag, and elongation lists).

Parameters:

cheby_fits (rubin_sim.movingObjects.chebyFits) – ChebyFits object, with attribute ‘coeffs’ - a dictionary of lists of coefficients.

class rubin_sim.moving_objects.DirectObs(footprint='camera', r_fov=1.75, x_tol=5, y_tol=3, eph_mode='nbody', prelim_eph_mode='2body', eph_type='basic', obs_code='I11', eph_file=None, obs_time_col='observationStartMJD', obs_time_scale='TAI', seeing_col='seeingFwhmGeom', visit_exp_time_col='visitExposureTime', obs_ra='fieldRA', obs_dec='fieldDec', obs_rot_sky_pos='rotSkyPos', obs_degrees=True, outfile_name='lsst_obs.dat', obs_info='', tstep=1.0, rough_tol=10.0, verbose=False, night_col='night', filter_col='filter', m5_col='fiveSigmaDepth', obs_id_col='observationId', pre_comp_tol=2.08)[source]

Bases: BaseObs

Generate observations of a set of moving objects: exact ephemeris at the times of each observation.

First generates observations on a rough grid and looks for observations within a specified tolerance of the actual observations; for the observations which pass this cut, generates a precise ephemeris and checks if the object is within the FOV.

Parameters:
  • footprint (str, optional) – Specify the footprint for the FOV. Options include “camera”, “circle”, “rectangle”. ‘Camera’ means use the actual LSST camera footprint (following a rough cut with a circular FOV). Default is circular FOV.

  • r_fov (float, optional) – If footprint is “circular”, this is the radius of the fov (in degrees). Default 1.75 degrees.

  • x_tol (float, optional) – If footprint is “rectangular”, this is half of the width of the (on-sky) fov in the RA direction (in degrees). Default 5 degrees.

  • y_tol (float, optional) – If footprint is “rectangular”, this is half of the width of the fov in Declination (in degrees). Default is 3 degrees

  • eph_mode (str, optional) – Mode for ephemeris generation - nbody or 2body. Default is nbody.

  • prelim_eph_mode (str, optional) – Mode for preliminary ephemeris generation, if any is done. Default is 2body.

  • eph_type (str, optional) – Type of ephemerides to generate - full or basic. Full includes all values calculated by openorb; Basic includes a more basic set. Default is Basic. (this includes enough information for most standard MAF metrics).

  • eph_file (str or None, optional) – The name of the planetary ephemerides file to use in ephemeris generation. Default (None) will use the default for PyOrbEphemerides.

  • obs_code (str, optional) – Observatory code for ephemeris generation. Default is “I11” - Cerro Pachon.

  • obs_time_col (str, optional) – Name of the time column in the obsData. Default ‘observationStartMJD’.

  • obs_time_scale (str, optional) – Type of timescale for MJD (TAI or UTC currently). Default TAI.

  • seeing_col (str, optional) – Name of the seeing column in the obsData. Default ‘seeingFwhmGeom’. This should be the geometric/physical seeing as it is used for the trailing loss calculation.

  • visit_exp_time_col (str, optional) – Name of the visit exposure time column in the obsData. Default ‘visitExposureTime’.

  • obs_ra (str, optional) – Name of the RA column in the obsData. Default ‘fieldRA’.

  • obs_dec (str, optional) – Name of the Dec column in the obsData. Default ‘fieldDec’.

  • obs_rot_sky_pos (str, optional) – Name of the Rotator column in the obsData. Default ‘rotSkyPos’.

  • obs_degrees (bool, optional) –

    Whether the observational data is in degrees or radians.

    Default True (degrees).

  • outfile_name (str, optional) – The output file name. Default is ‘lsst_obs.dat’.

  • obs_info (str, optional) – A string that captures provenance information about the observations. For example: ‘baseline_v2.0_10yrs, MJD 59853-61677’ or ‘baseline2018a minus NES’ Default ‘’.

  • tstep (float, optional) – The time between initial (rough) ephemeris generation points, in days. Default 1 day.

  • rough_tol (float, optional) – The initial rough tolerance value for positions, used as a first cut to identify potential observations (in degrees). Default 10 degrees.

  • pre_comp_tol (float (2.08)) – The radial tolerance to add when using pre-computed orbits. Should be larger than the full field of view extent.

run(orbits, obs_data, object_positions=None, object_mjds=None)[source]

Find and write the observations of each object to disk.

For each object, a rough grid of ephemeris points are either generated on the fly or read from a pre-calculated grid; If the rough grids indicate that an object may be present in an observation, then a more precise position is generated for the time of the observation.

Parameters:
  • orbits (rubin_sim.moving_objects.Orbits) – The orbits for which to generate ephemerides.

  • obs_data (np.ndarray) – The simulated pointing history data.

  • object_positions (np.ndarray) – Pre-computed RA,dec positions for each object in orbits (degrees).

  • object_mjds (np.ndarray) – MJD values for each pre-computed position.

class rubin_sim.moving_objects.Orbits[source]

Bases: object

Orbits reads, checks for required values, and stores orbit parameters for moving objects.

self.orbits stores the orbital parameters, as a pandas dataframe. self.dataCols defines the columns required, although obj_id, H, g, and sed_filename are optional.

assign_sed(orbits, random_seed=None)[source]

Assign either a C or S type SED, depending on the semi-major axis of the object. P(C type) = 0 (a<2); 0.5*a - 1 (2<a<4); 1 (a > 4), based on figure 23 from Ivezic et al 2001 (AJ, 122, 2749).

Parameters:

orbits (pd.DataFrame, pd.Series or np.ndarray) – Array-like object containing orbital parameter information.

Returns:

sedvals – Array containing the SED type for each object in ‘orbits’.

Return type:

np.ndarray

read_orbits(orbit_file, delim=None, skiprows=None)[source]

Read orbits from a file.

This generates a pandas dataframe containing columns matching dataCols, for the appropriate orbital parameter format. (currently accepts COM, KEP or CAR formats).

After reading and standardizing the column names, calls self.set_orbits to validate the orbital parameters. Expects angles in orbital element formats to be in degrees.

Note that readOrbits uses pandas.read_csv to read the data file with the orbital parameters. Thus, it should have column headers specifying the column names .. unless skiprows = -1 or there is just no header line at all. in which case it is assumed to be a standard DES format file, with no header line.

Parameters:
  • orbit_file (str) – The name of the input file with orbital parameter information.

  • delim (str, optional) – The delimiter for the input orbit file. Default is None, will use delim_whitespace=True.

  • skiprows (int, optional) – The number of rows to skip before reading the header information. Default is None, which will trigger a search of the file for the header columns.

set_orbits(orbits)[source]

Set and validate orbital parameters contain all required values.

Sets self.orbits and self.orb_format. If objid is not present in orbits, a sequential series of integers will be used. If H is not present in orbits, a default value of 20 will be used. If g is not present in orbits, a default value of 0.15 will be used. If sed_filename is not present in orbits, either C or S type will be assigned according to the semi-major axis value.

Parameters:

orbits (pd.DataFrame, pd.Series or np.ndarray) – Array-like object containing orbital parameter information.

update_orbits(neworb)[source]

Update existing orbits with new values, leaving OrbitIds, H, g, and sed_filenames in place.

Example use: transform orbital parameters (using PyOrbEphemerides) and then replace original values. Example use: propagate orbital parameters (using PyOrbEphemerides) and then replace original values.

Parameters:

neworb (pd.DataFrame)

class rubin_sim.moving_objects.PyOrbEphemerides(ephfile=None)[source]

Bases: object

Generate ephemerides and propagate orbits, using the python interface to Oorb.

PyOrbEphemerides handles the packing and unpacking of the fortran style arrays that pyoorb uses, to and from more user-friendly pandas arrays.

Parameters:

ephfile (str, optional) – Planetary ephemerides file for Oorb (i.e. de430 or de405). Default $OORB_DATA/de430.dat ($OORB_DATA = $OORB_DIR/data).

Examples

Typical usage:

>>> pyephs = PyOrbEphemerides()
>>> pyephs.set_orbits(orbits)
>>> ephs = pyephs.generateEphemerides(times, timeScale, obscode)
convert_from_oorb_elem()[source]

Translate pyoorb-style (fortran packed) orbital element array into a pandas dataframe. Operates on self.oorb_elem.

Returns:

new_orbits – A DataFrame with the appropriate subset of columns relating to orbital elements.

Return type:

pd.DataFrame

convert_orbit_format(orb_format='CAR')[source]

Convert orbital elements into format.

Example: converts from self.oorb_elem[orb_format] (such as KEP) to oorb_format (such as CAR).

Parameters:

orb_format (str, optional) – Format to convert orbital elements into.

generate_ephemerides(times, time_scale='UTC', obscode='I11', by_object=True, eph_mode='nbody', eph_type='basic')[source]

Calculate ephemerides for all orbits at times times.

All returned positions and angles are in degrees, velocities are degrees/day and distances are in AU.

Parameters:
  • times (np.ndarray, (N,)) – Ephemeris times.

  • time_scale (str, optional) – Time scale (UTC, TT, TAI) of times.

  • obscode (int or str, optional) – The observatory code for ephemeris generation.

  • by_object (bool, optional) – If True (default), resulting converted ephemerides are grouped by object. If False, resulting converted ephemerides are grouped by time.

  • eph_mode (str, optional) – Dynamical model to use for ephemeris generation - nbody or 2body. Accepts ‘nbody’, ‘2body’, ‘N’ or ‘2’. Default nbody.

  • eph_type (str, optional) – Generate full (more data) ephemerides or basic (less data) ephemerides. Default basic.

Returns:

ephemerides – The ephemeris values, organized as chosen by the user.

Return type:

np.ndarray

Notes

The returned ephemerides are a numpy array that can be grouped by object or by time.

If they are grouped by object (by_object = True), the array is organized as ephemeris_values[object][time]. Here the “ra” column is a 2-d array where the [0] axis length equals the number of ephemeris times.

If they are grouped by time (by_object=False), the array is organized as ephemeris_values[time][object]. Here the “ra” column is a 2-d array where the [0] axis length equals the number of objects.

propagate_orbits(new_epoch, eph_mode='nbody')[source]

Propagate orbits from self.orbits.epoch to new epoch (MJD TT).

Parameters:

new_epoch (float) – MJD TT time for new epoch.

set_orbits(orbit_obj)[source]

Set the orbits, to be used to generate ephemerides.

Immediately calls self._convertOorbElem to translate to the ‘packed’ oorb format.

Parameters:

orbit_obj (rubin_sim.moving_objects.Orbits) – The orbits to use to generate ephemerides.

rubin_sim.moving_objects.chebeval(x, p, interval=(-1.0, 1.0), do_velocity=True, mask=False)[source]

Evaluate a Chebyshev series and first derivative at points x.

If p is of length n + 1, this function returns: y_hat(x) = p_0 * T_0(x*) + p_1 * T_1(x*) + … + p_n * T_n(x*) where T_n(x*) are the orthogonal Chebyshev polynomials of the first kind, defined on the interval [-1, 1] and p_n are the coefficients. The scaled variable x* is defined on the [-1, 1] interval such that (x*) = (2*x - a - b)/(b - a), and x is defined on the [a, b] interval.

Parameters:
  • x (scalar or np.ndarray) – Points at which to evaluate the polynomial.

  • p (np.ndarray) – Chebyshev polynomial coefficients, as returned by chebfit.

  • interval (2-element list/tuple) – Bounds the x-interval on which the Chebyshev coefficients were fit.

  • do_velocity (bool) – If True, compute the first derivative at points x.

  • mask (bool) – If True, return Nans when the x goes beyond ‘interval’. If False, extrapolate fit beyond ‘interval’ limits.

Returns:

y, v – Y (position) and velocity values (if computed)

Return type:

float or np.ndarray, float or np.ndarray (or None)

rubin_sim.moving_objects.chebfit(t, x, dxdt=None, x_multiplier=None, dx_multiplier=None, n_poly=7)[source]

Fit Chebyshev polynomial constrained at endpoints using Newhall89 approach.

Return Chebyshev coefficients and statistics from fit to array of positions (x) and optional velocities (dx/dt). If both the function and its derivative are specified, then the value and derivative of the interpolating polynomial at the endpoints will be exactly equal to the input endpoint values. Many approximations may be piecewise strung together and the function value and its first derivative will be continuous across boundaries. If derivatives are not provided, only the function value will be continuous across boundaries.

If x_multiplier and dx_multiplier are not provided or are an inappropriate shape for t and x, they will be recomputed. See Newhall, X. X. 1989, Celestial Mechanics, 45, p. 305-310 for details.

Parameters:
  • t (np.ndarray) – Array of regularly sampled independent variable (e.g. time)

  • x (np.ndarray) – Array of regularly sampled dependent variable (e.g. declination)

  • dxdt (`np.ndarray’, optional) – Optionally, array of first derivatives of x with respect to t, at the same grid points. (e.g. sky velocity ddecl/dt)

  • x_multiplier (np.ndarray, optional) – Optional 2D Matrix with rows of C1^(-1)C2 corresponding to x. Use make_cheb_matrix to compute

  • dx_multiplier (np.ndarray, optional) – Optional 2D Matrix with rows of C1^(-1)C2 corresponding to dx/dt. Use make_cheb_matrix to compute

  • n_poly (int, optional) – Number of polynomial terms. Degree + 1. Must be >=2 and <=2*n_points, when derivative information is specified, or <=n_points, when no derivative information is specified. Default = 7.

Returns:

  • a_n (np.ndarray) – Array of chebyshev coefficients with length=n_poly.

  • residuals (np.ndarray) – Array of residuals of the tabulated function x minus the approximated function.

  • rms (float) – The rms of the residuals in the fit.

  • maxresid (float) – The maximum of the residals to the fit.

rubin_sim.moving_objects.get_oorb_data_dir()[source]

Find where the oorb files should be installed

rubin_sim.moving_objects.make_cheb_matrix(n_points, n_poly, weight=0.16)[source]

Compute C1^(-1)C2 using Newhall89 approach.

Utility function for fitting chebyshev polynomials to x(t) and dx/dt(t) forcing equality at the end points. This function computes the matrix (C1^(-1)C2). Multiplying this matrix by the x and dx/dt values to be fit produces the chebyshev coefficient. This function need only be called once for a given polynomial degree and number of points.

The matrices returned are of shape(n_points+1)x(n_poly). The coefficients fitting the n_points+1 points, X, are found by: A = xMultiplier * x + dxMultiplier * dxdt if derivative information is known, or A = xMultiplier * x if no derivative information is known. The xMultiplier matrices are different, depending on whether derivative information is known. Use function make_cheb_matrix_only_x if derviative is not known. See Newhall, X. X. 1989, Celestial Mechanics, 45, p. 305-310 for details.

Parameters:
  • n_points (int) – Number of point to be fits. Must be greater than 2.

  • n_poly (int) – Number of polynomial terms. Polynomial degree + 1

  • weight (float, optional) – Weight to allow control of relative effectos of position and velocity values. Newhall80 found best results are obtained with velocity weighted at 0.4 relative to position, giving W the form (1.0, 0.16, 1.0, 0.16,…)

Returns:

  • c1c2 (np.ndarray) – xMultiplier, C1^(-1)C2 even rows of shape (n_points+1)x(n_poly) to be multiplied by x values.

  • c1c2 (np.ndarray) – dxMultiplier, C1^(-1)C2 odd rows of shape (n_points+1)x(n_poly) to be multiplied by dx/dy values

rubin_sim.moving_objects.make_cheb_matrix_only_x(n_points, n_poly)[source]

Compute C1^(-1)C2 using Newhall89 approach without dx/dt

Compute xMultiplier using only the equality constraint of the x-values at the endpoints. To be used when first derivatives are not available. If chebyshev approximations are strung together piecewise only the x-values and not the first derivatives will be continuous at the boundaries. Multiplying this matrix by the x-values to be fit produces the chebyshev coefficients. This function need only be called once for a given polynomial degree and number of points. See Newhall, X. X. 1989, Celestial Mechanics, 45, p. 305-310.

Parameters:
  • n_points (int) – Number of point to be fits. Must be greater than 2.

  • n_poly (int) – Number of polynomial terms. Polynomial degree + 1

Returns:

c1c2 – xMultiplier, Even rows of C1^(-1)C2 w/ shape (n_points+1)x(n_poly) to be multiplied by x values

Return type:

np.ndarray

rubin_sim.moving_objects.read_observations(simfile, colmap, constraint=None, dbcols=None)[source]

Read the opsim database.

Parameters:
  • simfile (str) – Name (& path) of the opsim database file.

  • colmap (dict) – colmap dictionary (from rubin_sim.maf.batches.ColMapDict)

  • constraint (str, optional) – Optional SQL constraint (minus ‘where’) on the data to read from db. Default is None.

  • dbcols (list of [str], optional) – List of additional columns to query from the db and add to the output observations. Default None.

Returns:

simdata – The OpSim data read from the database.

Return type:

np.ndarray, (N)