lambertValidator#

Executive Summary#

This module computes the required Delta-V to change the spacecraft velocity to the one obtained by the lambertSolver module and checks if the resulting trajectory comes too close to the central body or too far away from the desired target location. The current spacecraft state from the NavTransMsgPayload navigation input message is propagated to the maneuver time \(t_{maneuver}\) using a 4th order Runge-Kutta (RK4) to obtain the expected state just before the Delta-V maneuver. The computed Delta-V and the specified state uncertainty are added to generate a set of initial states right after the maneuver. These initial states are propagated to the final time \(t_{final}\) where the spacecraft is supposed to arrive at the targeted position \({}^N\mathbf{r}_{T/N}\). If the final distance from the target position is less than the specified maximum distance \(r_{TB,max}\) for all trajectories, and none of the trajectories comes closer than \(r_{min}\) to the central body, the computed Delta-V is written to the DvBurnCmdMsgPayload output message. If any of these constraints are violated, the message content remains zero. The message is also zeroed if the LambertSolutionMsgPayload or LambertPerformanceMsgPayload messages indicate that the Lambert solution is not valid or has not converged yet.

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

navTransInMsg

NavTransMsgPayload

translational navigation input message

lambertProblemInMsg

lambertProblemMsgPayload

lambert problem setup input message

lambertSolutionInMsg

LambertSolutionMsgPayload

lambert problem solution input message

lambertPerformanceInMsg

LambertPerformanceMsgPayload

lambert problem performance message (additional information about the solution process)

dvBurnCmdOutMsg

DvBurnCmdMsgPayload

Delta-V command output message

Module Assumptions and Limitations#

The equations of motion used inside the module to propagate the state assume simple two body point mass gravity, and the motion is propagated using a 4th order Runge-Kutta (RK4). Additionally, this module assumes that \(t_{final} > t_{maneuver} > t\), with final time \(t_{final}\), maneuver time \(t_{maneuver}\) and current time \(t\).

Algorithm#

Equations of motion (two body point mass gravity) with gravitational parameter \(\mu\) and spacecraft position vector \(\mathbf{r}\):

(1)#\[\mathbf{\ddot{r}} = - \frac{\mu}{r^3} \mathbf{r}\]

The 27 perturbed initial states that are propagated to check for any constraint violations are obtained in the following way: The uncertainty of each state, specified by the 6x6 matrix “uncertaintyStates” in Hill frame components, is added to each state in the plus and minus direction. For \(N=6\) states, this gives \(2 \times N = 12\) initial states. The uncertainty of the Delta-V magnitude is applied to the Delta-V vector for each of those initial states, both in the Delta-V direction and against the Delta-V direction (corresponding to maximum and minimum expected DV magnitude). Applied to the 12 initial states, this gives 24 states. The last 3 initial states are obtained by applying the maximum and minimum expected DV to the unperturbed spacecraft state (+2), and by leaving the state and Delta-V vector entirely unperturbed (+1).

User Guide#

The module is first initialized as follows:

module = lambertValidator.LambertValidator()
module.modelTag = "lambertValidator"
module.finalTime = 2000.
module.maneuverTime = 1000.
module.maxDistanceTarget = 3000.
module.minOrbitRadius = 6378 * 1000.
module.uncertaintyStates = np.diag([5., 5., 5., 0.01, 0.01, 0.001])  # in Hill frame
module.uncertaintyDV = 0.1   # [m/s]
module.dvConvergenceTolerance = 0.01    # [m/s]
unitTestSim.AddModelToTask(unitTaskName, module)

The input messages are then connected:

module.navTransInMsg.subscribeTo(navTransInMsg)
module.lambertProblemInMsg.subscribeTo(lambertProblemInMsg)
module.lambertPerformanceInMsg.subscribeTo(lambertPerformanceInMsg)
module.lambertSolutionInMsg.subscribeTo(lambertSolutionInMsg)

Class LambertValidator#

class LambertValidator : public SysModel#

This module validates if the Lambert velocity solution violates any constraints and, if not, creates the Delta-V burn command message.

Public Functions

LambertValidator()#

This is the constructor for the module class. It sets default variable values and initializes the various parts of the model

~LambertValidator()#

Module Destructor

void reset(uint64_t currentSimNanos) override#

This method is used to reset the module and checks that required input messages are connected.

Parameters:

currentSimNanos – current simulation time in nano-seconds

Returns:

void

void updateState(uint64_t currentSimNanos) override#

This is the main method that gets called every time the module is updated.

Parameters:

currentSimNanos – current simulation time in nano-seconds

Returns:

void

Public Members

ReadFunctor<NavTransMsgPayload> navTransInMsg#

translational navigation input message

ReadFunctor<LambertProblemMsgPayload> lambertProblemInMsg#

lambert problem input message

ReadFunctor<LambertSolutionMsgPayload> lambertSolutionInMsg#

lambert solution input message

ReadFunctor<LambertPerformanceMsgPayload> lambertPerformanceInMsg#

lambert performance input message

Message<DvBurnCmdMsgPayload> dvBurnCmdOutMsg#

Delta-V burn command message.

Message<LambertValidatorMsgPayload> lambertValidatorOutMsg#

Lambert Validator results message.

BSKLogger bskLogger#

&#8212; BSK Logging

double lambertSolutionSpecifier = 1#

[-] which Lambert solution (1 or 2), if applicable, should be used

double finalTime = {}#

[s] time at which target position should be reached

double maneuverTime = {}#

[s] time at which maneuver should be executed

double maxDistanceTarget = {}#

[m] maximum acceptable distance from target location at final time

double minOrbitRadius = {}#

6x6 matrix square root of the covariance matrix to apply errors with, in Hill (Orbit) frame components

[m] minimum acceptable orbit radius

double uncertaintyDV = 0.1#

[m/s] uncertainty of the Delta-V magnitude

double dvConvergenceTolerance = 1e-3#

[m/s] tolerance on difference between DeltaV solutions between time steps

bool ignoreConstraintViolations = false#

override flag to write DV message despite constraint violations