Module: orbElemOffset
Executive Summary
The OrbElemOffset module combines two classical orbital element messages into a single output element set. One
input provides the nominal (main) orbit, and the second provides element offsets. By default, all fields are added
element wise, including true anomaly. Optionally, the anomaly offset can be interpreted as a mean anomaly increment.
Module Description
The module reads:
mainElementsInMsg: baseline classical elementsoffsetElementsInMsg: element offsets to apply
For the following fields, the output is always computed through direct addition:
The handling of f depends on useMeanAnomalyOffset:
False(default): treatoffset.fas a true-anomaly increment\[f_{out} = f_{main} + f_{off}\]True: treatoffset.fas a mean-anomaly increment \(\Delta M\)Convert \((f_{main}, e_{main})\) to \(M_{main}\)
Form \(e_{out} = e_{main} + e_{off}\)
Form \(M_{out} = M_{main} + \Delta M\)
Convert \((M_{out}, e_{out})\) back to \(f_{out}\)
This preserves a mean-anomaly interpretation for the offset while still allowing eccentricity to be updated through
offset.e.
Message Interfaces
Msg Variable Name |
Msg Type |
Description |
|---|---|---|
mainElementsInMsg |
Input message containing the baseline classical elements
( |
|
offsetElementsInMsg |
Input message containing additive offsets for the classical elements. The |
|
elementsOutMsg |
Output message containing the combined classical elements after applying offsets. |
Detailed Behavior
At each update step, the module performs the following operations:
Reads baseline elements from
mainElementsInMsg.Reads element offsets from
offsetElementsInMsg.Computes
a,e,i,Omega, andomegavia element wise addition.Computes
feither as a direct true-anomaly addition or via mean-anomaly conversion/update based onuseMeanAnomalyOffset.Writes the resulting
ClassicElementsMsgPayloadtoelementsOutMsg.
During Reset(), the module checks that both input messages are linked and logs an error if either input is not
connected.
Verification and Testing
The module is verified by an automated unit test in
src/fswAlgorithms/orbitControl/orbElemOffset/_UnitTest/test_orbElemOffset.py.
The test is parameterized over both anomaly handling modes:
useMeanAnomalyOffset = False: verifies all six fields are direct sums.useMeanAnomalyOffset = True: verifies the non-anomaly fields are direct sums and verifies that the output mean anomaly differs from the baseline mean anomaly by exactlyoffset.fwithin numerical tolerance.
This validates both the baseline additive behavior and the mean-anomaly-offset conversion path.
-
class OrbElemOffset : public SysModel
- #include <orbElemOffset.h>
Classic orbital elements combiner with optional mean anomaly offset handling.
This module reads two ClassicElementsMsgPayload inputs:
mainElementsInMsg: nominal or “base” classical orbital elements
offsetElementsInMsg: element offsets to apply to the base elements
For most fields the module performs a simple element wise addition
out.a = main.a + offset.a out.e = main.e + offset.e out.i = main.i + offset.i out.Omega = main.Omega + offset.Omega out.omega = main.omega + offset.omega
The treatment of the true anomaly field f depends on the flag useMeanAnomalyOffset
If useMeanAnomalyOffset is false (default) then the offsetElementsInMsg.f field is interpreted as a true anomaly increment and the output is
out.f = main.f + offset.f
If useMeanAnomalyOffset is true then offsetElementsInMsg.f is interpreted as a mean anomaly offset \(\Delta M\) in radians The module then
1 Converts the main elements (main.f main.e) to mean anomaly M_main 2 Forms the final eccentricity as e_out = main.e + offset.e 3 Forms the final mean anomaly M_out = M_main + offset.f 4 Converts (M_out e_out) back to a true anomaly f_out
and sets out.f = f_out
All other elements in ClassicElementsMsgPayload are left to zero.
The output payload is written to elementsOutMsg.
Public Functions
-
OrbElemOffset() = default
Default constructor
-
void Reset(uint64_t CurrentSimNanos) override
Reset the internal state and validate input message connections.
This method checks that both mainElementsInMsg and offsetElementsInMsg are connected. If either is not linked an error is logged through the logger.
- Parameters:
CurrentSimNanos – Current simulation time in nanoseconds (unused)
-
void UpdateState(uint64_t CurrentSimNanos) override
Compute the combined classical orbital elements and write them out.
The module reads both input messages, constructs the output classical elements according to the rules described in the class documentation, and writes the result to elementsOutMsg
- Parameters:
CurrentSimNanos – Current simulation time in nanoseconds
Public Members
-
bool useMeanAnomalyOffset = false
Toggle interpretation of offsetElementsInMsg.f as a mean anomaly offset.
false: offsetElementsInMsg.f is treated as a true anomaly increment
true: offsetElementsInMsg.f is treated as a mean anomaly offset \(\Delta M\)
-
ReadFunctor<ClassicElementsMsgPayload> mainElementsInMsg
Nominal or “main” classical elements input message.
-
ReadFunctor<ClassicElementsMsgPayload> offsetElementsInMsg
Classical elements offset input message.
-
Message<ClassicElementsMsgPayload> elementsOutMsg
Output classical elements message for the combined elements.
-
BSKLogger bskLogger
BSK Logging.