Module: orbitalElementControl
Executive Summary
The OrbitalElementControl module specializes LinearTimeInvariantSystem to
implement orbital-element feedback control that outputs an inertial force command. It reads target and current
classical orbital elements, forms a six-element error vector, applies an LTI law, maps the result through the
classical Gauss control matrix, and publishes CmdForceInertialMsgPayload.
Module Description
The inherited LTI dynamics are
with fixed input and output size for this subclass:
input size: \(m = 6\)
output size: \(p = 6\)
The input vector is constructed from target and current classical elements as
where:
\(\Delta a / a = (a_{current} - a_{target})/a_{target}\)
\(\Delta e = e_{current} - e_{target}\)
\(\Delta i, \Delta\Omega, \Delta\omega\) are wrapped into \([-\pi,\pi]\)
\(\Delta M\) is formed by converting each true anomaly to mean anomaly and wrapping the difference into \([-\pi,\pi]\)
The LTI output \(\mathbf{y}\) is interpreted in orbital-element space and mapped to LVLH force using the classical Gauss matrix \(\mathbf{B}_{oe}\):
Then LVLH force is rotated into inertial coordinates and written to
forceOutMsg.forceRequestInertial.
Convenience gain helpers are provided:
setProportionalGain(K): assigns \(\mathbf{D}=K\) (6x6)setIntegralGain(K): assigns \(\mathbf{B}=I_{6}\) and \(\mathbf{C}=K\) (6x6), enabling integral action through the internal LTI state
Message Interfaces
Msg Variable Name |
Msg Type |
Description |
|---|---|---|
targetOEInMsg |
Input target classical orbital elements used as the controller reference. |
|
currentOEInMsg |
Input current classical orbital elements used to form feedback errors and to build the Gauss mapping and frame transforms. |
|
forceOutMsg |
Output inertial-frame commanded force vector in |
Detailed Behavior
At each update step, the module performs the following operations:
Reads
targetOEInMsgandcurrentOEInMsg.Constructs the six-element input error vector \(\mathbf{u}\) with wrapped angular terms and mean-anomaly difference.
Uses inherited LTI dynamics to compute \(\dot{\mathbf{x}}\) and \(\mathbf{y}\).
Builds the classical Gauss matrix \(\mathbf{B}_{oe}\) from current orbital elements and configured
mu.Computes LVLH force \(\mathbf{f}_{LVLH} = -\mathbf{B}_{oe}^{T}\mathbf{y}\).
Converts current orbital elements to inertial position/velocity, builds the LVLH-to-inertial rotation, and rotates \(\mathbf{f}_{LVLH}\) into inertial force.
Publishes inertial force through
forceOutMsg.
During Reset(), the module:
runs base-class matrix-size consistency checks,
verifies
mu > 0,verifies both
targetOEInMsgandcurrentOEInMsgare linked.
Verification and Testing
The module is verified in
src/simulation/mujocoDynamics/linearTimeInvariantSystem/_UnitTest/test_orbitalElementControl.py.
The test is parameterized for both proportional and integral gain configurations and compares
OrbitalElementControl output against meanOEFeedback configured with effectively zero \(J_2\). For the same
initial chief/deputy orbit setup and gains, the inertial force command from OrbitalElementControl is required to
match the reference meanOEFeedback force within tolerance.
-
class OrbitalElementControl : public LinearTimeInvariantSystem
- #include <orbitalElementControl.h>
LTI based orbital element feedback controller that outputs an inertial force command.
This class derives from LinearTimeInvariantSystem and implements a controller in classical orbital element coordinates. The LTI input is a 6 element vector of orbital element errors
u = [da/a, de, di, dOmega, domega, dM]^T
constructed from target and current classical elements. The LTI output y is also a 6 element vector in the same orbital element space.
The module then uses the Gauss planetary control matrix B(oe) to map y to a Hill frame force, and finally transforms this force into the inertial frame using the current orbit geometry. The inertial force is written to a CmdForceInertialMsgPayload.
The internal state dimension and A, B, C, D matrices are configured through the LinearTimeInvariantSystem interface and are not fixed by this class.
Public Functions
-
OrbitalElementControl() = default
Default constructor.
Constructs an OrbitalElementControl instance with mu set to zero and unconnected input and output messages. The user must set mu, configure the LTI matrices, and connect the input and output messages before the module is used in a simulation.
-
void setProportionalGain(const Eigen::MatrixXd &K)
Set the proportional gain matrix for the orbital element controller.
This method configures the proportional feedback gain used in the internal LinearTimeInvariantSystem realization of the control law.
The proportional term maps the instantaneous 6x1 orbital-element error vector into a 6x1 control output before it is transformed by the Gauss planetary-equations B-matrix in writeOutput().
Internally, this method assigns the gain matrix K to the D matrix of the underlying LTI system:
y = C x + D u
where:
u is the six-element orbital-element error vector
D = K implements the proportional action
- Parameters:
K – 6x6 proportional gain matrix.
-
void setIntegralGain(const Eigen::MatrixXd &K)
Set the integral gain matrix for the orbital element controller.
This method configures the integral feedback gain for the controller, implemented using the state dynamics of the underlying LinearTimeInvariantSystem. The integral term accumulates the orbital-element error over time and contributes to the output force.
The structure set here is:
ẋ = B u with B = I [6x6] y = C x + D u with C = K_I
This realizes:
x(t) = ∫ u(t) dt (integrator state) y = K_I x (integral feedback)
- Parameters:
K – 6x6 integral gain matrix.
-
virtual void Reset(uint64_t CurrentSimNanos) override
Reset the module and validate configuration.
This method first calls LinearTimeInvariantSystem::Reset to validate A, B, C, and D against the sizes returned by getStateSize, getInputSize, and getOutputSize. It then checks that the gravitational parameter mu is positive and that both targetOEInMsg and currentOEInMsg are linked.
- Parameters:
CurrentSimNanos – Current simulation time in nanoseconds.
-
inline virtual size_t getInputSize() const override
Size of the LTI input vector.
The input to the LTI system is the 6 element orbital element error vector
u = [da/a, de, di, dOmega, domega, dM]^T.
- Returns:
Number of input elements, always 6.
-
inline virtual size_t getOutputSize() const override
Size of the LTI output vector.
The LTI output y is a 6 element vector in orbital element space. It is not directly the applied force. Instead, writeOutput uses y together with the Gauss control matrix B to compute a Hill frame force, then transforms that force into the inertial frame and writes it to the forceOutMsg.
- Returns:
Number of output elements, always 6.
-
virtual Eigen::VectorXd readInput(uint64_t CurrentSimNanos) override
Read the current input vector for the LTI system.
This method reads the target and current classical orbital element messages and constructs the 6 element orbital element error vector
u = [da/a, de, di, dOmega, domega, dM]^T,
where:
da/a is the relative semi major axis error
de is the eccentricity error
di is the inclination error
dOmega is the RAAN error
domega is the argument of periapsis error
dM is the mean anomaly error
Angular differences are wrapped into [-pi, pi].
- Parameters:
CurrentSimNanos – Current simulation time in nanoseconds.
- Returns:
Input vector u as an Eigen::VectorXd of length getInputSize().
-
virtual void writeOutput(uint64_t CurrentSimNanos, const Eigen::VectorXd &y) override
Map the LTI output vector to an inertial force command and write the output message.
The argument y is the output of the LTI system in orbital element space, with dimension getOutputSize() = 6. This method performs the following steps:
Read the current classical orbital elements from currentOEInMsg and construct a ClassicElements struct.
Compute the 6 by 3 Gauss control matrix B(oe) that maps Hill frame force to orbital element rates.
Interpret y as an orbital elements feedback vector and compute a Hill frame force, for example via f_H = B^T * y.
Convert the current elements to inertial position and velocity and compute the Hill to inertial direction cosine matrix.
Rotate the Hill frame force into the inertial frame and write it into forceOutMsg as a CmdForceInertialMsgPayload.
- Parameters:
CurrentSimNanos – Current simulation time in nanoseconds.
y – LTI output vector in orbital element space with length getOutputSize().
Public Members
-
double mu = 0.0
Gravitational parameter mu [m^3/s^2].
This parameter must be set to the central body’s gravitational parameter before use. The Reset method checks that mu is positive and reports an error if it is not.
-
ReadFunctor<ClassicElementsMsgPayload> targetOEInMsg
Target classical orbital elements input.
-
ReadFunctor<ClassicElementsMsgPayload> currentOEInMsg
Current classical orbital elements input.
-
Message<CmdForceInertialMsgPayload> forceOutMsg
Inertial force command output message.
-
OrbitalElementControl() = default