grss.libgrss.PropSimulation#

class grss.libgrss.PropSimulation#

Bases: pybind11_object

The PropSimulation class contains the orbit propagation simulation for intgrating solar system small bodies.

__init__(*args, **kwargs)#

Overloaded function.

  1. __init__(self: grss.libgrss.PropSimulation, name: str, t0: float, defaultSpiceBodies: int, DEkernelPath: str) -> None

    Constructor for the PropSimulation class.

    namestr

    Name of the simulation.

    t0real

    Initial MJD TDB time of the simulation.

    defaultSpiceBodiesint

    Version of the DE kernel to get the default SPICE bodies from.

    DEkernelPathstr

    Path to the SPICE DE kernel.

  2. __init__(self: grss.libgrss.PropSimulation, name: str, simRef: grss.libgrss.PropSimulation) -> None

    Constructor for the PropSimulation class.

    namestr

    Name of the simulation.

    simRefPropSimulation

    Simulation to copy.

Methods

__init__(*args, **kwargs)

Overloaded function.

add_event(self, event)

Adds an event to the simulation.

add_integ_body(self, body)

Adds an integration body to the simulation.

add_spice_body(self, body)

Adds a SPICE body to the simulation.

extend(self, tf[, tEvalNew, observerInfoNew])

Extends the simulation to a new final time.

get_integration_parameters(self)

Gets the integration parameters.

get_sim_constants(self)

Gets the constants of the simulation.

get_spiceBody_state(self, t, bodyName)

Gets the state of a SPICE body at a given time.

integrate(self)

Propagates the simulation using the Gauss-Radau integrator.

interpolate(self, t)

Interpolates the states of the integrated bodies to a given time.

map_ephemeris(self)

Memory maps the ephemeris of the simulation.

remove_body(self, name)

Removes a body from the simulation.

save(self, filename)

Saves the simulation to a file.

set_integration_parameters(self, tf[, ...])

Sets the integration parameters.

set_sim_constants(self[, du2m, tu2s, G, clight])

Sets the constants of the simulation.

unmap_ephemeris(self)

Unmaps the ephemeris of the simulation.

Attributes

DEkernelPath

Path to the SPICE DE kernel.

caParams

Close approach parameters of the simulation.

consts

Constants of the simulation.

convergedLightTime

Whether to use converged Newtonian light time correction.

evalApparentState

Whether to evaluate the apparent state of the integration bodies.

evalMeasurements

Whether to evaluate the measurements of the integration bodies.

impactParams

Impact parameters of the simulation.

integBodies

Integration bodies of the simulation.

integParams

Integration parameters of the simulation.

interpParams

Interpolation parameters of the simulation.

lightTimeEval

Light time from the observer to each integration body for each value in PropSimulation.tEval.

name

Name of the simulation.

obsType

Observation type for each value in PropSimulation.tEval (0=optical, 1=delay, 2=doppler, 3=Gaia).

observerInfo

Observer information for each value in PropSimulation.tEval.

opticalObs

Optical observation of each integration body in the simulation for each value in PropSimulation.tEval.

opticalObsCorr

Photocenter-barycenter correction for each optical observation for each integration body in the simulation.

opticalObsDot

Time derivative of the optical observation of each integration body in the simulation for each value in PropSimulation.tEval.

opticalPartials

Optical observation partials of each integration body in the simulation for each value in PropSimulation.tEval.

radarObs

Radar observation of each integration body in the simulation for each value in PropSimulation.tEval.

radarPartials

Radar observation partials of each integration body in the simulation for each value in PropSimulation.tEval.

spiceBodies

SPICE bodies of the simulation.

t

Current time of the simulation.

tEval

MJD Times to evaluate the states of the integrated bodies at.

tEvalMargin

Margin for allowing evaluation past the propagation start and end times.

tEvalUTC

Whether the MJD evaluation time is in UTC for each value in PropSimulation.tEval, as opposed to TDB.

unsafePersistentMemoryMap

Whether to use unsafe persistent memory mapping for the simulation.

xInteg

Current states of each integration body in the simulation.

xIntegEval

States of each integration body in the simulation for each value in PropSimulation.tEval.

xObserver

State of the observer for each value in PropSimulation.tEval.

property DEkernelPath#

Path to the SPICE DE kernel.

add_event(self: grss.libgrss.PropSimulation, event: grss.libgrss.Event) None#

Adds an event to the simulation.

Parameters:

event (libgrss.Event) – Event to add to the simulation.

add_integ_body(self: grss.libgrss.PropSimulation, body: grss.libgrss.IntegBody) None#

Adds an integration body to the simulation.

Parameters:

body (libgrss.IntegBody) – Integration body to add to the simulation.

add_spice_body(self: grss.libgrss.PropSimulation, body: grss.libgrss.SpiceBody) None#

Adds a SPICE body to the simulation.

Parameters:

body (libgrss.SpiceBody) – SPICE body to add to the simulation.

property caParams#

Close approach parameters of the simulation. List of libgrss.CloseApproachParameters objects.

property consts#

Constants of the simulation. libgrss.Constants object.

property convergedLightTime#

Whether to use converged Newtonian light time correction.

property evalApparentState#

Whether to evaluate the apparent state of the integration bodies.

property evalMeasurements#

Whether to evaluate the measurements of the integration bodies.

extend(self: grss.libgrss.PropSimulation, tf: float, tEvalNew: list[float] = [], observerInfoNew: list[list[float]] = []) None#

Extends the simulation to a new final time.

Parameters:
  • tf (real) – New final time of integration (MJD TDB).

  • tEvalNew (list of real) – Extra MJD Times to evaluate the states of the integrated bodies at. Can be TDB or UTC based on tEvalUTC.

  • observerInfoNew (list of list of real) – New observer information. Each list at least contains the central body SPICE ID (e.g., 399 for Earth) and the body-fixed longitude, latitude, and distance. This information should be repeated for radar observations.

get_integration_parameters(self: grss.libgrss.PropSimulation) list[float]#

Gets the integration parameters.

Returns:

  • nInteg (int) – Number of integrated bodies.

  • nSpice (int) – Number of bodies with SPICE ephemerides.

  • nTotal (int) – Total number of bodies. nTotal = nInteg + nSpice.

  • t0 (real) – Initial time of integration (MJD TDB).

  • tf (real) – Final time of integration (MJD TDB).

  • adaptiveTimestep (bool) – Flag to use adaptive time step for the propagation.

  • dt0 (real) – Initial time step.

  • dtMin (real) – Minimum time step.

  • dtChangeFactor (real) – Factor by which to limit the change in time step.

  • tolInteg (real) – Tolerance for integration.

  • tolPC (real) – Tolerance for predictor-corrector within IAS15.

get_sim_constants(self: grss.libgrss.PropSimulation) list[float]#

Gets the constants of the simulation.

Returns:

  • du2m (real) – Conversion factor from distance units to meters.

  • tu2s (real) – Conversion factor from time units to seconds.

  • G (real) – Gravitational constant.

  • clight (real) – Speed of light in a vacuum.

  • j2000Jd (real) – Julian date of J2000 epoch.

  • JdMinusMjd (real) – Difference between Julian date and modified Julian date.

get_spiceBody_state(self: grss.libgrss.PropSimulation, t: float, bodyName: str) list[float]#

Gets the state of a SPICE body at a given time.

Parameters:
  • t (real) – Time to get the state at.

  • bodyName (str) – Name of the SPICE body in the simulation.

property impactParams#

Impact parameters of the simulation. List of libgrss.ImpactParameters objects.

property integBodies#

Integration bodies of the simulation. List of libgrss.IntegBody objects.

property integParams#

Integration parameters of the simulation. libgrss.IntegParams object.

integrate(self: grss.libgrss.PropSimulation) None#

Propagates the simulation using the Gauss-Radau integrator.

property interpParams#

Interpolation parameters of the simulation. libgrss.InterpolationParameters object.

interpolate(self: grss.libgrss.PropSimulation, t: float) list[float]#

Interpolates the states of the integrated bodies to a given time.

Parameters:

t (real) – Time to interpolate to.

Returns:

xIntegInterp – Interpolated GEOMETRIC states of the integrated bodies.

Return type:

list of real

property lightTimeEval#

Light time from the observer to each integration body for each value in PropSimulation.tEval.

map_ephemeris(self: grss.libgrss.PropSimulation) None#

Memory maps the ephemeris of the simulation.

Returns:

None – None.

Return type:

NoneType

property name#

Name of the simulation.

property obsType#

Observation type for each value in PropSimulation.tEval (0=optical, 1=delay, 2=doppler, 3=Gaia).

property observerInfo#

Observer information for each value in PropSimulation.tEval.

property opticalObs#

Optical observation of each integration body in the simulation for each value in PropSimulation.tEval.

property opticalObsCorr#

Photocenter-barycenter correction for each optical observation for each integration body in the simulation.

property opticalObsDot#

Time derivative of the optical observation of each integration body in the simulation for each value in PropSimulation.tEval.

property opticalPartials#

Optical observation partials of each integration body in the simulation for each value in PropSimulation.tEval.

property radarObs#

Radar observation of each integration body in the simulation for each value in PropSimulation.tEval.

property radarPartials#

Radar observation partials of each integration body in the simulation for each value in PropSimulation.tEval.

remove_body(self: grss.libgrss.PropSimulation, name: str) None#

Removes a body from the simulation.

Parameters:

name (str) – Name of the body to remove.

save(self: grss.libgrss.PropSimulation, filename: str) None#

Saves the simulation to a file.

Parameters:

filename (str) – Name of the file to save the simulation to.

set_integration_parameters(self: grss.libgrss.PropSimulation, tf: float, tEval: list[float] = [], tEvalUTC: bool = False, evalApparentState: bool = False, convergedLightTime: bool = False, observerInfo: list[list[float]] = [], adaptiveTimestep: bool = True, dt0: float = 1.0, dtMin: float = 0.0001, dtChangeFactor: float = 0.25, tolInteg: float = 1e-11, tolPC: float = 1e-16) None#

Sets the integration parameters.

Parameters:
  • tf (real) – Final time of integration (MJD TDB).

  • tEval (list of real) – MJD Times to evaluate the states of the integrated bodies at. Can be TDB or UTC based on tEvalUTC.

  • tEvalUTC (bool) – Whether the MJD evaluation time is in UTC for each value in libgrss.tEval, as opposed to TDB.

  • evalApparentState (bool) – Whether to evaluate the apparent state of the integration bodies.

  • convergedLightTimes (bool) – Whether to use converged Newtonian light time correction.

  • observerInfo (list of list of real) – Observer information. Each list at least contains the central body SPICE ID (e.g., 399 for Earth) and the body-fixed longitude, latitude, and distance. This information should be repeated for radar observations.

  • adaptiveTimestep (bool) – Flag to use adaptive time step for the propagation.

  • dt0 (real) – Initial time step.

  • dtMin (real) – Minimum time step.

  • dtChangeFactor (real) – Factor by which to limit the change in time step.

  • tolInteg (real) – Tolerance for integration.

  • tolPC (real) – Tolerance for predictor-corrector within IAS15.

set_sim_constants(self: grss.libgrss.PropSimulation, du2m: float = 149597870700.0, tu2s: float = 86400.0, G: float = 1.4881851702345195e-34, clight: float = 173.14463267424034) None#

Sets the constants of the simulation.

Parameters:
  • du2m (real) – Conversion factor from distance units to meters.

  • tu2s (real) – Conversion factor from time units to seconds.

  • G (real) – Gravitational constant.

  • clight (real) – Speed of light in a vacuum.

property spiceBodies#

SPICE bodies of the simulation. List of libgrss.SpiceBodies objects.

property t#

Current time of the simulation.

property tEval#

MJD Times to evaluate the states of the integrated bodies at. Can be TDB or UTC based on PropSimulation.tEvalUTC.

property tEvalMargin#

Margin for allowing evaluation past the propagation start and end times.

property tEvalUTC#

Whether the MJD evaluation time is in UTC for each value in PropSimulation.tEval, as opposed to TDB.

unmap_ephemeris(self: grss.libgrss.PropSimulation) None#

Unmaps the ephemeris of the simulation.

Returns:

None – None.

Return type:

NoneType

property unsafePersistentMemoryMap#

Whether to use unsafe persistent memory mapping for the simulation.

property xInteg#

Current states of each integration body in the simulation.

property xIntegEval#

States of each integration body in the simulation for each value in PropSimulation.tEval.

property xObserver#

State of the observer for each value in PropSimulation.tEval.