lsst.sims.movingObjects package

Submodules

lsst.sims.movingObjects.baseObs module

class lsst.sims.movingObjects.baseObs.BaseObs(footprint='circle', rFov=1.75, xTol=5, yTol=3, ephMode='nbody', ephType='basic', obsCode='I11', ephFile=None, obsTimeCol='observationStartMJD', obsTimeScale='TAI', seeingCol='seeingFwhmGeom', visitExpTimeCol='visitExposureTime', obsRA='fieldRA', obsDec='fieldDec', obsRotSkyPos='rotSkyPos', obsDegrees=True, outfileName='lsst_obs.dat', obsMetadata='')[source]

Bases: object

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

Parameters:
  • footPrint (str, opt) – 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.
  • rFov (float, opt) – If footprint is “circular”, this is the radius of the fov (in degrees). Default 1.75 degrees.
  • xTol (float, opt) – If footprint is “rectangular”, this is half of the width of the (on-sky) fov in the RA direction (in degrees). Default 5 degrees. (so size of footprint in degrees will be 10 degrees in the RA direction).
  • yTol (float, opt) – If footprint is “rectangular”, this is half of the width of the fov in Declination (in degrees). Default is 3 degrees (so size of footprint in degrees will be 6 degrees in the Dec direction).
  • ephMode (str, opt) – Mode for ephemeris generation - nbody or 2body. Default is nbody.
  • ephType (str, opt) – 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).
  • ephFile (str or None, opt) – The name of the planetary ephemerides file to use for ephemeris generation. Default (None) will use the default for PyOrbEphemerides.
  • obsCode (str, opt) – Observatory code for ephemeris generation. Default is “I11” - Cerro Pachon.
  • obsTimeCol (str, opt) – Name of the time column in the obsData. Default ‘observationStartMJD’.
  • obsTimeScale (str, opt) – Type of timescale for MJD (TAI or UTC currently). Default TAI.
  • seeingCol (str, opt) – 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.
  • visitExpTimeCol (str, opt) – Name of the visit exposure time column in the obsData. Default ‘visitExposureTime’.
  • obsRA (str, opt) – Name of the RA column in the obsData. Default ‘fieldRA’.
  • obsDec (str, opt) – Name of the Dec column in the obsData. Default ‘fieldDec’.
  • obsRotSkyPos (str, opt) – Name of the Rotator column in the obsData. Default ‘rotSkyPos’.
  • obsDegrees (bool, opt) – Whether the observational data is in degrees or radians. Default True (degrees).
  • outfileName (str, opt) – The output file name. Default is ‘lsst_obs.dat’.
  • obsMetadata (str, opt) – A string that captures provenance information about the observations. For example: ‘kraken_2026, MJD 59853-61677’ or ‘baseline2018a minus NES’ Default ‘’.
calcColors(sedname='C.dat', sedDir=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 (opt)) – Name of the SED. Default ‘C.dat’.
  • sedDir (str (opt)) – Directory containing the SEDs of the moving objects. Default None = $SIMS_MOVINGOBJECTS_DIR/data.
Returns:

Dictionary of the colors in self.filterlist for this particular Sed.

Return type:

dict

calcTrailingLosses(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). dmagDetect 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, opt) – The exposure time of the images, in seconds. Default 30.
Returns:

dmagTrail and dmagDetect for each set of velocity/seeing/texp values.

Return type:

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

generateEphemerides(sso, times, ephMode=None, ephType=None)[source]

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

The default engine here is pyoorb, however this method could be overwritten to use another ephemeris generator, such as ADAM.

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

Parameters:
  • sso (lsst.sims.movingObjects.Orbits) – Typically this will be a single object.
  • times (np.ndarray) – The times at which to generate ephemerides. MJD.
  • ephMode (str or None, opt) – Potentially override default ephMode (self.ephMode). Must be ‘2body’ or ‘nbody’.
Returns:

Ephemerides of the sso.

Return type:

pandas.Dataframe

readFilters(filterDir=None, bandpassRoot='total_', bandpassSuffix='.dat', filterlist=('u', 'g', 'r', 'i', 'z', 'y'), vDir=None, vFilter='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:
  • filterDir (str, opt) – Directory containing the filter throughput curves (‘total*.dat’) Default set by ‘LSST_THROUGHPUTS_BASELINE’ env variable.
  • bandpassRoot (str, opt) – Rootname of the throughput curves in filterlist. E.g. throughput curve names are bandpassRoot + filterlist[i] + bandpassSuffix Default ‘total_’ (appropriate for LSST throughput repo).
  • bandpassSuffix (str, opt) – Suffix for the throughput curves in filterlist. Default ‘.dat’ (appropriate for LSST throughput repo).
  • filterlist (list, opt) – List containing the filter names to use to calculate colors. Default (‘u’, ‘g’, ‘r’, ‘i’, ‘z’, ‘y’)
  • vDir (str, opt) – Directory containing the V band throughput curve. Default None = $SIMS_MOVINGOBJECTS_DIR/data.
  • vFilter (str, opt) – Name of the V band filter curve. Default harris_V.dat.
setupEphemerides()[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.

ssoInCameraFov(ephems, obsData)[source]

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

Parameters:
  • ephems (np.recarray) – Ephemerides for the objects.
  • obsData (np.recarray) – Observation pointings.
Returns:

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

Return type:

np.ndarray

ssoInCircleFov(ephems, obsData)[source]

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

Parameters:
  • ephems (np.recarray) – Ephemerides for the objects.
  • obsData (np.recarray) – The observation pointings.
Returns:

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

Return type:

numpy.ndarray

ssoInFov(ephems, obsData)[source]

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

Parameters:
  • ephems (np.recarray) – Ephemerides for the objects.
  • obsData (np.recarray) – Observation pointings.
Returns:

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

Return type:

np.ndarray

ssoInRectangleFov(ephems, obsData)[source]

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

Parameters:
  • ephems (np.recarray) – Ephemerides for the objects.
  • obsData (np.recarray) – The observation pointings.
Returns:

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

Return type:

numpy.ndarray

writeObs(objId, objEphs, obsData, sedname='C.dat')[source]

Call for each object; write out the observations of each object.

This method is called once all of the ephemeris values for each observation are calculated; the calling method should have already done the matching between ephemeris & simulated observations to find the observations where the object is within the specified fov. Inside this method, the trailing losses and color terms are calculated and added to the output observation file.

The first time this method is called, a header will be added to the output file.

Parameters:
  • objId (str or int or float) – The identifier for the object (from the orbital parameters)
  • objEphs (numpy.ndarray) – The ephemeris values of the object at each observation. Note that the names of the columns are encoded in the numpy structured array, and any columns included in the returned ephemeris array will also be propagated to the output.
  • obsData (numpy.ndarray) – The observation details from the simulated pointing history, for all observations of the object. All columns automatically propagated to the output file.
  • sedname (str, out) – The sed_filename for the object (from the orbital parameters). Used to calculate the appropriate color terms for the output file. Default “C.dat”.

lsst.sims.movingObjects.chebyFits module

class lsst.sims.movingObjects.chebyFits.ChebyFits(orbitsObj, tStart, tSpan, timeScale='TAI', obscode=807, skyTolerance=2.5, nCoeff_position=14, nCoeff_vmag=9, nCoeff_delta=5, nCoeff_elongation=6, ngran=64, ephFile=None, nDecimal=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. 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 nDecimal 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 for non-database use; 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 skyTolerance - 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.

Parameters:
  • orbitsObj (Orbits) – The orbits for which to fit chebyshev polynomial coefficients.
  • tStart (float) – The starting point in time to fit coefficients. MJD.
  • tSpan (float) – The time span (starting at tStart) over which to fit coefficients. Days.
  • timeScale ({'TAI', 'UTC', 'TT'}) – The timescale of the MJD time, tStart, and the timeScale that should be used with the chebyshev coefficients.
  • obsCode (int, optional) – The observatory code of the location for which to generate ephemerides. Default 807 (CTIO).
  • skyTolerance (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.
  • ephFile (str, optional) – The path to the JPL ephemeris file to use. Default is ‘$OORB_DATA/de405.dat’.
  • nDecimal (int, optional) – The number of decimal places to allow in the segment length (and thus the times of the endpoints) can be limited to nDecimal places. Default 10. For LSST SIMS moving object database, this should be 13 decimal places for NEOs and 0 for all others.
calcOneSegment(orbitObj, 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:
  • orbitObj (Orbits) – The single Orbits object we’re fitting at the moment.
  • ephs (numpy.ndarray) – The ephemerides we’re fitting at the moment (for the single object / single segment).
calcSegmentLength(length=None)[source]

Set the typical initial ephemeris timestep and segment length for all objects between tStart/tEnd.

Sets self.length.

The segment length will fit into the time period between tStart/tEnd 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.
calcSegments()[source]

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

generateEphemerides(times, byObject=True)[source]

Generate ephemerides using OpenOrb for all orbits.

Parameters:times (numpy.ndarray) – The times to use for ephemeris generation.
makeAllTimes()[source]

Using tStart and tEnd, 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:Numpy array of times.
Return type:numpy.ndarray
write(coeffFile, residFile, failedFile, append=False)[source]

Write coefficients, residuals and failed fits to disk.

Parameters:
  • coeffFile (str) – The filename for the coefficient values.
  • residFile (str) – The filename for the residual values.
  • failedFile (str) – The filename to write the failed fit information (if failed objects exist).
  • append (bool, optional) – Flag to append (or overwrite) the output files.

lsst.sims.movingObjects.chebyValues module

class lsst.sims.movingObjects.chebyValues.ChebyValues[source]

Bases: object

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

getEphemerides(times, objIds=None, extrapolate=False)[source]

Find the ephemeris information for ‘objIds’ 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.
  • objIds (numpy.ndarray, optional) – The object ids for which to generate ephemerides. If None, then just uses all objects.
  • extrapolate (bool) – 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:

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

Return type:

numpy.ndarray

readCoefficients(chebyFitsFile)[source]

Read coefficients from output file written by ChebyFits.

Parameters:chebyFitsFile (str) – The filename of the coefficients file.
setCoefficients(chebyFits)[source]

Set coefficients using a ChebyFits object. (which contains a dictionary of objId, tStart, tEnd, ra, dec, delta, vmag, and elongation lists).

Parameters:chebyFits (chebyFits) – ChebyFits object, with attribute ‘coeffs’ - a dictionary of lists of coefficients.

lsst.sims.movingObjects.chebyshevUtils module

Utilities to estimate and evaluate Chebyshev coefficients of a function.

Implementation of Newhall, X. X. 1989, Celestial Mechanics, 45, p. 305-310

lsst.sims.movingObjects.chebyshevUtils.chebeval(x, p, interval=(-1.0, 1.0), doVelocity=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 numpy.ndarray) – Points at which to evaluate the polynomial.
  • p (numpy.ndarray) – Chebyshev polynomial coefficients, as returned by chebfit.
  • interval (2-element list/tuple) – Bounds the x-interval on which the Chebyshev coefficients were fit.
  • doVelocity (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 (position) and velocity values (if computed)

Return type:

scalar or numpy.ndarray, scalar or numpy.ndarray

lsst.sims.movingObjects.chebyshevUtils.chebfit(t, x, dxdt=None, xMultiplier=None, dxMultiplier=None, nPoly=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 xMultiplier and dxMultiplier 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 (numpy.ndarray) – Array of regularly sampled independent variable (e.g. time)
  • x (numpy.ndarray) – Array of regularly sampled dependent variable (e.g. declination)
  • dxdt (numpy.ndarray, optional) – Optionally, array of first derivatives of x with respect to t, at the same grid points. (e.g. sky velocity ddecl/dt)
  • xMultiplier (numpy.ndarray, optional) – Optional 2D Matrix with rows of C1^(-1)C2 corresponding to x. Use makeChebMatrix to compute
  • dxMultiplier (numpy.ndarray, optional) – Optional 2D Matrix with rows of C1^(-1)C2 corresponding to dx/dt. Use makeChebMatrix to compute
  • nPoly (int, optional) – Number of polynomial terms. Degree + 1. Must be >=2 and <=2*nPoints, when derivative information is specified, or <=nPoints, when no derivative information is specified. Default = 7.
Returns:

  • numpy.ndarray – Array of chebyshev coefficients with length=nPoly.
  • numpy.ndarray – Array of residuals of the tabulated function x minus the approximated function.
  • float – The rms of the residuals in the fit.
  • float – The maximum of the residals to the fit.

lsst.sims.movingObjects.chebyshevUtils.makeChebMatrix(nPoints, nPoly, 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(nPoints+1)x(nPoly). The coefficients fitting the nPoints+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 makeChebMatrixOnlyX if derviative is not known. See Newhall, X. X. 1989, Celestial Mechanics, 45, p. 305-310 for details.

Parameters:
  • nPoints (int) – Number of point to be fits. Must be greater than 2.
  • nPoly (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:

  • numpy.ndarray – xMultiplier, C1^(-1)C2 even rows of shape (nPoints+1)x(nPoly) to be multiplied by x values.
  • numpy.ndarray – dxMultiplier, C1^(-1)C2 odd rows of shape (nPoints+1)x(nPoly) to be multiplied by dx/dy values

lsst.sims.movingObjects.chebyshevUtils.makeChebMatrixOnlyX(nPoints, nPoly)[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:
  • nPoints (int) – Number of point to be fits. Must be greater than 2.
  • nPoly (int) – Number of polynomial terms. Polynomial degree + 1
Returns:

xMultiplier, Even rows of C1^(-1)C2 w/ shape (nPoints+1)x(nPoly) to be multiplied by x values

Return type:

numpy.ndarray

lsst.sims.movingObjects.directObs module

class lsst.sims.movingObjects.directObs.DirectObs(footprint='circular', rFov=1.75, xTol=5, yTol=3, ephMode='nbody', prelimEphMode='2body', ephType='basic', obsCode='I11', ephFile=None, obsTimeCol='observationStartMJD', obsTimeScale='TAI', seeingCol='seeingFwhmGeom', visitExpTimeCol='visitExposureTime', obsRA='fieldRA', obsDec='fieldDec', obsRotSkyPos='rotSkyPos', obsDegrees=True, outfileName='lsst_obs.dat', obsMetadata='', tstep=1.0, roughTol=10.0)[source]

Bases: lsst.sims.movingObjects.baseObs.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, opt) – 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.
  • rFov (float, opt) – If footprint is “circular”, this is the radius of the fov (in degrees). Default 1.75 degrees.
  • xTol (float, opt) – If footprint is “rectangular”, this is half of the width of the (on-sky) fov in the RA direction (in degrees). Default 5 degrees. (so size of footprint in degrees will be 10 degrees in the RA direction).
  • yTol (float, opt) – If footprint is “rectangular”, this is half of the width of the fov in Declination (in degrees). Default is 3 degrees (so size of footprint in degrees will be 6 degrees in the Dec direction).
  • ephMode (str, opt) – Mode for ephemeris generation - nbody or 2body. Default is nbody.
  • prelimEphMode (str, opt) – Mode for preliminary ephemeris generation, if any is done. Default is 2body.
  • ephType (str, opt) – 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).
  • ephFile (str or None, opt) – The name of the planetary ephemerides file to use for ephemeris generation. Default (None) will use the default for PyOrbEphemerides.
  • obsCode (str, opt) – Observatory code for ephemeris generation. Default is “I11” - Cerro Pachon.
  • obsTimeCol (str, opt) – Name of the time column in the obsData. Default ‘observationStartMJD’.
  • obsTimeScale (str, opt) – Type of timescale for MJD (TAI or UTC currently). Default TAI.
  • seeingCol (str, opt) – 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.
  • visitExpTimeCol (str, opt) – Name of the visit exposure time column in the obsData. Default ‘visitExposureTime’.
  • obsRA (str, opt) – Name of the RA column in the obsData. Default ‘fieldRA’.
  • obsDec (str, opt) – Name of the Dec column in the obsData. Default ‘fieldDec’.
  • obsRotSkyPos (str, opt) – Name of the Rotator column in the obsData. Default ‘rotSkyPos’.
  • obsDegrees (bool, opt) – Whether the observational data is in degrees or radians. Default True (degrees).
  • outfileName (str, opt) – The output file name. Default is ‘lsst_obs.dat’.
  • obsMetadata (str, opt) – A string that captures provenance information about the observations. For example: ‘kraken_2026, MJD 59853-61677’ or ‘baseline2018a minus NES’ Default ‘’.
  • tstep (float, opt) – The time between initial (rough) ephemeris generation points, in days. Default 1 day.
  • roughTol (float, opt) – The initial rough tolerance value for positions, used as a first cut to identify potential observations (in degrees). Default 10 degrees.
run(orbits, obsData)[source]

Find and write the observations of each object to disk.

For each object, generate a very rough grid of ephemeris points (typically using 2body integration). Then identify pointings in obsData which are within

Parameters:
  • orbits (lsst.sims.movingObjects.Orbits) – The orbits to generate ephemerides for.
  • obsData (np.recarray) – The simulated pointing history data.

lsst.sims.movingObjects.linearObs module

class lsst.sims.movingObjects.linearObs.LinearObs(footprint='circular', rFov=1.75, xTol=5, yTol=3, ephMode='nbody', ephType='basic', obsCode='I11', ephFile=None, obsTimeCol='observationStartMJD', obsTimeScale='TAI', seeingCol='seeingFwhmGeom', visitExpTimeCol='visitExposureTime', obsRA='fieldRA', obsDec='fieldDec', obsRotSkyPos='rotSkyPos', obsDegrees=True, outfileName='lsst_obs.dat', obsMetadata='', tstep=0.08333333333333333)[source]

Bases: lsst.sims.movingObjects.baseObs.BaseObs

Generate observations for a set of Orbits using linear interpolation.

Uses linear interpolations between grid of true ephemerides. Ephemerides can be generated using 2-body or n-body integration.

Parameters:
  • footPrint (str, opt) – 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.
  • rFov (float, opt) – If footprint is “circle”, this is the radius of the fov (in degrees). Default 1.75 degrees.
  • xTol (float, opt) – If footprint is “rectangle”, this is half of the width of the (on-sky) fov in the RA direction (in degrees). Default 5 degrees. (so size of footprint in degrees will be 10 degrees in the RA direction).
  • yTol (float, opt) – If footprint is “rectangular”, this is half of the width of the fov in Declination (in degrees). Default is 3 degrees (so size of footprint in degrees will be 6 degrees in the Dec direction).
  • ephMode (str, opt) – Mode for ephemeris generation - nbody or 2body. Default is nbody.
  • ephType (str, opt) – 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).
  • ephFile (str or None, opt) – The name of the planetary ephemerides file to use for ephemeris generation. Default (None) will use the default for PyOrbEphemerides.
  • obsCode (str, opt) – Observatory code for ephemeris generation. Default is “I11” - Cerro Pachon.
  • obsTimeCol (str, opt) – Name of the time column in the obsData. Default ‘observationStartMJD’.
  • obsTimeScale (str, opt) – Type of timescale for MJD (TAI or UTC currently). Default TAI.
  • seeingCol (str, opt) – 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.
  • visitExpTimeCol (str, opt) – Name of the visit exposure time column in the obsData. Default ‘visitExposureTime’.
  • obsRA (str, opt) – Name of the RA column in the obsData. Default ‘fieldRA’.
  • obsDec (str, opt) – Name of the Dec column in the obsData. Default ‘fieldDec’.
  • obsRotSkyPos (str, opt) – Name of the Rotator column in the obsData. Default ‘rotSkyPos’.
  • obsDegrees (bool, opt) – Whether the observational data is in degrees or radians. Default True (degrees).
  • outfileName (str, opt) – The output file name. Default is ‘lsst_obs.dat’.
  • obsMetadata (str, opt) – A string that captures provenance information about the observations. For example: ‘kraken_2026, MJD 59853-61677’ or ‘baseline2018a minus NES’ Default ‘’.
  • tstep (float, opt) – The time between points in the ephemeris grid, in days. Default 2 hours.
run(orbits, obsData)[source]

Find and write the observations of each object to disk.

For each object, identify the observations where the object is within rFOV of the pointing boresight (potentially, also in the camera footprint), and write the ephemeris values and observation metadata to disk. Uses linear interpolation between ephemeris gridpoints.

Parameters:
  • orbits (lsst.sims.movingObjects.Orbits) – The orbits to generate ephemerides for.
  • obsData (np.recarray) – The simulated pointing history data.

lsst.sims.movingObjects.lsstCameraFootprint module

class lsst.sims.movingObjects.lsstCameraFootprint.LsstCameraFootprint[source]

Bases: object

Class to provide the capability for identifying observations within an LSST camera footprint.

inCameraFov(ephems, obsData, epoch=2000.0, timeCol='observationStartMJD')[source]

Determine which observations are within the actual camera footprint for a series of observations.

Parameters:
  • ephems (np.recarray) – Ephemerides for the objects, with RA and Dec as ‘ra’ and ‘dec’ columns (in degrees).
  • obsData (np.recarray) – Observation pointings, with RA and Dec as ‘ra’ and ‘dec’ columns (in degrees). The telescope rotation angle should be in ‘rotSkyPos’ (in degrees), and the time of each pointing should be in the ‘expMJD’ column.
  • epoch (float, opt) – The epoch of the ephemerides and pointing data. Default 2000.0.
Returns:

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

Return type:

np.ndarray

lsst.sims.movingObjects.ooephemerides module

class lsst.sims.movingObjects.ooephemerides.PyOrbEphemerides(ephfile=None)[source]

Bases: object

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

Typical usage: pyephs = PyOrbEphemerides() # Set the orbital parameters, using an lsst.sims.movingObjects.Orbits object pyephs.setOrbits(orbits) # Generate ephemerides at times ‘times’. ephs = pyephs.generateEphemerides(times, timeScale=’UTC’, obscode=’I11’)

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

Parameters:ephfile (str, opt) – Planetary ephemerides file for Oorb (i.e. de430 or de405). Default $OORB_DATA/de430.dat ($OORB_DATA = $OORB_DIR/data).
convertFromOorbElem()[source]

Translate pyoorb-style orbital element array back into dataframe.

Parameters:oorbElem (numpy.ndarray) – The orbital elements in OpenOrb format.
Returns:A DataFrame with the appropriate subset of columns relating to orbital elements.
Return type:pd.DataFrame
convertOrbitFormat(orb_format='CAR')[source]

Convert orbital elements from the format in orbitObj into ‘format’.

Parameters:format (str, opt) – Format to convert orbital elements into.
generateEphemerides(times, timeScale='UTC', obscode='I11', byObject=True, ephMode='nbody', ephType='basic')[source]

Calculate ephemerides for all orbits at times times.

This is a public method, wrapping self._convertTimes, self._generateOorbEphs and self._convertOorbEphs (which include dealing with oorb-formatting of arrays).

The return ephemerides are in a numpy recarray, with axes - if byObject = True : [ephemeris values][object][@time] (i.e. the ‘ra’ column = 2-d array, where the [0] axis (length) equals the number of ephTimes) - if byObject = False : [ephemeris values][time][@object] (i.e. the ‘ra’ column = 2-d arrays, where the [0] axis (length) equals the number of objects)

The ephemeris values returned to the user (== columns of the recarray) are: [‘delta’, ‘ra’, ‘dec’, ‘magV’, ‘time’, ‘dradt’, ‘ddecdt’, ‘phase’, ‘solarelon’, ‘velocity’] where positions/angles are all in degrees, velocities are deg/day, and delta is the distance between the Earth and the object in AU.

Parameters:
  • ephtimes (numpy.ndarray) – Ephemeris times in oorb format (see self.convertTimes)
  • obscode (int or str, optional) – The observatory code for ephemeris generation. Default=807 (Cerro Tololo).
  • byObject (boolean, optional) – If True (default), resulting converted ephemerides are grouped by object. If False, resulting converted ephemerides are grouped by time.
  • ephMode (str, optional) – Dynamical model to use for ephemeris generation - nbody or 2body. Accepts ‘nbody’, ‘2body’, ‘N’ or ‘2’. Default nbody.
  • ephType (str, optional) – Generate full (more data) ephemerides or basic (less data) ephemerides. Default basic.
Returns:

The ephemeris values, organized as chosen by the user.

Return type:

numpy.ndarray

propagateOrbits(newEpoch)[source]

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

Parameters:new_epoch (float) – MJD TT time for new epoch.
setOrbits(orbitObj)[source]

Set the orbits, to be used to generate ephemerides.

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

Parameters:orbitObj (Orbits) – The orbits to use to generate ephemerides.

lsst.sims.movingObjects.orbits module

class lsst.sims.movingObjects.orbits.Orbits[source]

Bases: object

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

Instantiate the class and then use readOrbits or setOrbits to set the orbit values.

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

assignSed(orbits, randomSeed=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 (pandas.DataFrame, pandas.Series or numpy.ndarray) – Array-like object containing orbital parameter information.
Returns:Array containing the SED type for each object in ‘orbits’.
Return type:numpy.ndarray
readOrbits(orbitfile, delim=None, skiprows=None)[source]

Read orbits from a file, generating 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.setOrbits 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:
  • orbitfile (str) – The name of the input file containing 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 for pandas. Default is None, which will trigger a check of the file to look for the header columns.
setOrbits(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 (pandas.DataFrame, pandas.Series or numpy.ndarray) – Array-like object containing orbital parameter information.
updateOrbits(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 (pandas.DataFrame) –

lsst.sims.movingObjects.version module

Module contents

Moving object utilities.