ephemDifferenceWithUncertainty#

Executive Summary#

Module reads in two ephemeris messages containing the position and velocity of the two corresponding objects (Base and Secondary). The module then computes the relative position and velocity of the Secondary with respect to the Base, and writes this to a translational navigation message. Additionally, the module takes the 6x6 state covariance of the two bodies as input parameters via setter functions, and writes the relative state vector as well as the relative covariance matrix to a filter message.

The time tag that is written to the navigation and filter messages is taken from the Secondary (the base ephemeris may not be updated as frequently).

Message Connection Descriptions#

The following table lists all the module input and output messages. The module msg connection is set by the user from python. The msg type contains a link to the message structure definition, while the description provides information on what this message is used for:

Module I/O Messages#

Msg Variable Name

Msg Type

Description

ephemBaseInMsg

EphemerisMsgPayload

Input ephemeris message containing the position and velocity of the base object in the inertial frame

ephemSecondaryInMsg

EphemerisMsgPayload

Input ephemeris message containing the position and velocity of the secondary object in the inertial frame

navTransOutMsg

NavTransMsgPayload

Output navigation message containing the relative position and velocity in the inertial frame

filterOutMsg

FilterMsgPayload

Output filter message containing the relative state and covariance matrix

Detailed Module Description#

The relative position and velocity of the Secondary with respect to the Base are computed by

\[\begin{split}\mathbf{r}_{21}^N &= \mathbf{r}_{2}^N - \mathbf{r}_{1}^N \\ \mathbf{v}_{21}^N &= \mathbf{v}_{2}^N - \mathbf{v}_{1}^N\end{split}\]

where \(\mathbf{r}_{1}^N\) and \(\mathbf{v}_{1}^N\) are the position and velocity vectors of the base object, and \(\mathbf{r}_{2}^N\) and \(\mathbf{v}_{2}^N\) are the position and velocity vectors of the secondary object, respectively. The relative state vector is

\[\mathbf{X}_{21} = \mathbf{X}_{2} - \mathbf{X}_{1} = [(\mathbf{r}_{21}^N)^T, (\mathbf{v}_{21}^N)^T]^T\]

The relative covariance is found by:

\[\begin{split}Cov(\mathbf{X}_{21}, \mathbf{X}_{21}) &= Cov(\mathbf{X}_{2}, \mathbf{X}_{2}) - Cov(\mathbf{X}_{2}, \mathbf{X}_{1}) - Cov(\mathbf{X}_{1}, \mathbf{X}_{2}) + Cov(\mathbf{X}_{1}, \mathbf{X}_{1}) \\ &= Cov(\mathbf{X}_{1}, \mathbf{X}_{1}) + Cov(\mathbf{X}_{2}, \mathbf{X}_{2})\end{split}\]

where \(Cov(\mathbf{X}_{1}, \mathbf{X}_{1})\) and \(Cov(\mathbf{X}_{2}, \mathbf{X}_{2})\) are the covariance matrices of the base and secondary, respectively. Because the two objects are assumed to be independent, \(Cov(\mathbf{X}_{2}, \mathbf{X}_{1}) = Cov(\mathbf{X}_{1}, \mathbf{X}_{2}) = 0\).

User Guide#

This section is to outline the steps needed to setup the module in Python.

  1. Import the module:

    from Basilisk.fswAlgorithms import ephemDifferenceWithUncertainty
    
  2. Create an instantiation of the class:

    module = ephemDifferenceWithUncertainty.EphemDifferenceWithUncertainty()
    
  3. The covariance matrices are set by:

    module.setCovarianceBase(covar_1)
    module.setCovarianceSecondary(covar_2)
    
  4. Subscribe to the messages:

    module.ephemBaseInMsg.subscribeTo(ephem1InMsg)
    module.ephemSecondaryInMsg.subscribeTo(ephem2InMsg)
    
  5. Add model to task:

    sim.AddModelToTask(taskName, module)
    

Class EphemDifferenceWithUncertainty#

class EphemDifferenceWithUncertainty : public SysModel#

This module computes the difference between two ephemeris messages, and outputs the relative states into a navigation and filter message.

Public Functions

void updateState(uint64_t currentSimNanos) override#

During an update, this module computes the difference between two ephemeris messages.

Parameters:

currentSimNanos – The clock time at which the function was called (nanoseconds)

Returns:

void

void reset(uint64_t currentSimNanos) override#

This method resets the module state to its initialized default.

Parameters:

currentSimNanos – The clock time at which the function was called (nanoseconds)

Returns:

void

void setCovarianceBase(const Eigen::MatrixXd stateCovariance)#

Set the state covariance of the base celestial object (e.g. asteroid)

Parameters:

Eigen::MatrixXd – covariance

Returns:

void

Eigen::MatrixXd getCovarianceBase() const#

Get the state covariance of the base celestial object (e.g. asteroid)

Returns:

Eigen::MatrixXd covariance

void setCovarianceSecondary(const Eigen::MatrixXd stateCovariance)#

Set the state covariance of the secondary celestial object (e.g. spacecraft)

Parameters:

Eigen::MatrixXd – covariance

Returns:

void

Eigen::MatrixXd getCovarianceSecondary() const#

Get the state covariance of the secondary celestial object (e.g. spacecraft)

Returns:

Eigen::MatrixXd covariance