Module: meanOEFeedback
Summary
The primary purpose of this module is to calculate feedback control input for deputy spacecraft in terms of mean orbital element difference. This module uses Lyapunov control theory described in chapter 14 of Analytical Mechanics of Space Systems. In addition to classic orbital element set \((a,e,i,\Omega,\omega,M)\), this module can also deal with equinoctial orbital element \((a,e\sin{(\omega+\Omega)},e\cos{(\omega+\Omega)},\tan{(i/2)}\sin{(\Omega)},\tan{(i/2)}\cos{(\Omega)},\Omega+\omega+M)\) in order to avoid singularity. A control input \(\bf{u}\) used in this module is described in Eq. (1).
Control matrix \(B(oe)\) derived from Gauss’s Planetary Equation is given in Eq.(14.213) of Analytical Mechanics of Space Systems for classic orbital element, or in Application of several control techniques for the ionospheric observation nanosatellite formation for equinoctial one. Depending on preference, you can change which orbital element set to use by setting \(oe_{type}\) parameter. In order to calculate mean oe, position and velocity input messages are converted to osculating oe. Then, osculating oe is converted to mean oe by Brower’s theory.
Message Connection Descriptions
The following table lists all the module input and output messages. The module msg variable name 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.
Msg Variable Name |
Msg Type |
Description |
---|---|---|
chiefTransInMsg |
The name of the chief’s position and velocity input message |
|
deputyTransInMsg |
The name of the deputy’s position and velocity input message |
|
forceOutMsg |
Calculated Force to control orbital element difference output message |
Module Assumptions and Limitations
This module assumes that target orbital element differnce is constant during simulation.
User Guide
This module requires the following variables to be set as parameters:
oeType
0 for classic oe (default), 1 for equinoctial oemu
gravitational constant for a central body in m^3/s^2req
equatorial radius of a central body in metersJ2
J2 constant of a central bodytargetDiffOeMean
desired mean orbital element difference. Assumed oe is linked tooeType
.K
control gain matrix
You must be careful about a meter unit used for mu
and req
.
For targetDiffOeMean
, normalized semi major axis must be used.
Functions
-
void SelfInit_meanOEFeedback(meanOEFeedbackConfig *configData, int64_t moduleID)
This method initializes the configData for this module. It checks to ensure that the inputs are sane and then creates the output message
- Parameters:
configData – The configuration data associated with this module
moduleID – The Basilisk module identifier
- Returns:
void
-
void Update_meanOEFeedback(meanOEFeedbackConfig *configData, uint64_t callTime, int64_t moduleID)
Add a description of what this main Update() routine does for this module
- Parameters:
configData – The configuration data associated with the module
callTime – The clock time at which the function was called (nanoseconds)
moduleID – The Basilisk module identifier
- Returns:
void
-
void Reset_meanOEFeedback(meanOEFeedbackConfig *configData, uint64_t callTime, int64_t moduleID)
This method performs a complete reset of the module. Local module variables that retain time varying states between function calls are reset to their default values. The local copy of the message output buffer should be cleared.
- Parameters:
configData – The configuration data associated with the module
callTime – The clock time at which the function was called (nanoseconds)
moduleID – The Basilisk module identifier
- Returns:
void
-
struct meanOEFeedbackConfig
- #include <meanOEFeedback.h>
Top level structure for the sub-module routines.
Public Members
-
NavTransMsg_C chiefTransInMsg
chief orbit input message
-
NavTransMsg_C deputyTransInMsg
deputy orbit input message
-
CmdForceInertialMsg_C forceOutMsg
deputy control force output message
-
double K[36]
Lyapunov Gain (6*6)
-
double targetDiffOeMean[6]
target mean orbital element difference
-
uint8_t oeType
0: classic (default), 1: equinoctial
-
double mu
[m^3/s^2] gravitational constant
-
double req
[m] equatorial planet radius
-
double J2
[] J2 planet oblateness parameter
-
BSKLogger *bskLogger
BSK Logging.
-
NavTransMsg_C chiefTransInMsg