Module: gravityEffector¶
Executive Summary¶
Abstract class that is used to implement an effector impacting a GRAVITY body that does not itself maintain a state or represent a changing component of the body (for example: gravity, thrusters, solar radiation pressure, etc.)
The module
PDF Description
contains further information on this module’s function,
how to run it, as well as testing.
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.
Msg Variable Name |
Msg Type |
Description |
---|---|---|
centralBodyOutMsg |
central planet body state output message |
The gravity effector contains a list of GravBodyData
objects which contain the planet mass and size properties etc.
The following table lists the Spice planet ephemeris input message that can be connected to a GravBodyData
object.
If no message is connected, then the planet has zero position and orientation information by default.
Msg Variable Name |
Msg Type |
Description |
---|---|---|
planetBodyInMsg |
planet spice ephemerisis input message |
User Guide¶
The user must provide a list of GravBodyData
objects to the spacecraft using:
scObject.gravField.gravBodies = spacecraft.GravBodyVector(gravBodyList)
Each gravity body data object can be created using:
earth = gravityEffector.GravBodyData()
earth.planetName = 'earth_planet_data'
earth.mu = 0.3986004415E+15 # meters^3/s^2
earth.radEquator = 6378136.6 # meters
earth.isCentralBody = False
earth.useSphericalHarmParams = False
Note that the simIncludeGradBody.py
helper file contains a gravity body factor class to facilitate
setting up gravity bodies.
-
class Polyhedral¶
- #include <gravityEffector.h>
poyhedral class
Public Functions
-
Polyhedral()¶
-
~Polyhedral()¶
-
bool initializeParameters()¶
[-] configure polyhedral based on inputs
-
Eigen::Vector3d computeField(const Eigen::Vector3d pos_Pfix)¶
———————————Main Interface————————-—///
Use to compute the field in position pos, given in a body frame.
- Parameters
pos_Pfix – Position in which the field is to be computed.
- Returns
acc Vector including the computed field.
-
bool polyReady()¶
class variable
Public Members
-
unsigned int nVertex¶
[-] Number of vertexes
-
unsigned int nFacet¶
[-] Number of facets
-
double volPoly¶
[-] Volume of the polyhedral
-
double muBody¶
[-] Gravitation parameter for the planet
-
Eigen::MatrixXd xyzVertex¶
[m] Position of vertex
-
Eigen::MatrixXd orderFacet¶
[-] Vertexes of a facet
-
Eigen::MatrixXd normalFacet¶
[-] Normal of a facet
-
BSKLogger bskLogger¶
— BSK Logging
-
Polyhedral()¶
-
class SphericalHarmonics¶
- #include <gravityEffector.h>
spherical harmonics class
Public Functions
-
SphericalHarmonics()¶
-
~SphericalHarmonics()¶
-
bool initializeParameters()¶
[-] configure all spher-harm based on inputs
-
double getK(const unsigned int degree)¶
class method
-
Eigen::Vector3d computeField(const Eigen::Vector3d pos_Pfix, unsigned int degree, bool include_zero_degree)¶
———————————Main Interface————————-—///
Use to compute the field in position pos, given in a body frame.
- Parameters
pos_Pfix – Position in which the field is to be computed.
degree – used to compute the field.
include_zero_degree – Boolean that determines whether the zero-degree term is included.
- Returns
acc Vector including the computed field.
-
bool harmReady()¶
class variable
Public Members
-
unsigned int maxDeg¶
[-] Maximum degree of the spherical harmonics
-
double radEquator¶
[-] Reference radius for the planet
-
double muBody¶
[-] Gravitation parameter for the planet
-
std::vector<std::vector<double>> cBar¶
[-] C coefficient set
-
std::vector<std::vector<double>> sBar¶
[-] S coefficient set
-
std::vector<std::vector<double>> aBar¶
[-] Normalized ‘derived’ Assoc. Legendre
-
std::vector<std::vector<double>> n1¶
[-] What am I
-
std::vector<std::vector<double>> n2¶
[-] What am I
-
std::vector<std::vector<double>> nQuot1¶
[-] What am I
-
std::vector<std::vector<double>> nQuot2¶
[-] What am I
-
BSKLogger bskLogger¶
— BSK Logging
-
SphericalHarmonics()¶
-
class GravBodyData¶
- #include <gravityEffector.h>
Container for gravitational body data.
This class is designed to hold all of the information for a gravity body. The nominal use-case has it initialized at the python level and attached to dynamics using the AddGravityBody method.
Public Functions
-
GravBodyData()¶
Set parameters for a gravity body.
-
~GravBodyData()¶
Destructor.
-
void initBody(int64_t moduleID)¶
Method to initialize the gravity body.
-
Eigen::Vector3d computeGravityInertial(Eigen::Vector3d r_I, uint64_t simTimeNanos)¶
compute the gravitational acceleration
- Parameters
r_I – inertial position vector
simTimeNanos – simulation time (ns)
-
double computePotentialEnergy(Eigen::Vector3d r_I)¶
compute potential energy
- Parameters
r_I – inertial position vector
-
void loadEphemeris()¶
Command to load the ephemeris data.
load ephemeris information
- Parameters
moduleID –
-
void registerProperties(DynParamManager &statesIn)¶
class method
Public Members
-
ReadFunctor<SpicePlanetStateMsgPayload> planetBodyInMsg¶
planet spice ephemeris input message
-
bool isCentralBody = 0¶
Flag indicating that object is center.
-
bool usePolyhedral = 0¶
Flag indicating to use polyhedral model.
-
bool useSphericalHarmParams = 0¶
Flag indicating to use spherical harmonics perturbations.
-
double mu = 0¶
[m3/s^2] central body gravitational param
-
double ephemTime¶
[s] Ephemeris time for the body in question
-
double ephIntTime¶
[s] Integration time associated with the ephem data
-
double radEquator = 0¶
[m] Equatorial radius for the body
-
double radiusRatio = 0¶
[] ratio of polar over equatorial radius
-
SpicePlanetStateMsgPayload localPlanet = {}¶
[-] Class storage of ephemeris info from scheduled portion
-
uint64_t timeWritten = 0¶
[ns] time the input planet state message was written
-
std::string planetName = ""¶
Gravitational body name, this is used as the Spice name if spiceInterface is used.
-
std::string displayName = ""¶
this is the name that is displayed in Vizard. If not set, Vizard shows planetName
-
std::string modelDictionaryKey = ""¶
“” will result in using the current default for the celestial body’s given name, otherwise key will be matched if possible to available model in internal model dictionary
-
Polyhedral poly¶
Object that computes the polyhedral gravity field.
-
SphericalHarmonics spherHarm¶
Object that computes the spherical harmonics gravity field.
-
BSKLogger bskLogger¶
— BSK Logging
-
Eigen::MatrixXd *r_PN_N¶
[m] (state engine property) planet inertial position vector
-
Eigen::MatrixXd *v_PN_N¶
[m/s] (state engine property) planet inertial velocity vector
-
Eigen::MatrixXd *muPlanet¶
[m/s] (state engine property) planet inertial velocity vector
-
Eigen::MatrixXd *J20002Pfix¶
[m/s] (state engine property) planet attitude [PN]
-
Eigen::MatrixXd *J20002Pfix_dot¶
[m/s] (state engine property) planet attitude rate [PN_dot]
-
GravBodyData()¶
-
class GravityEffector : public SysModel¶
- #include <gravityEffector.h>
gravity effector class
Public Functions
-
GravityEffector()¶
-
~GravityEffector()¶
-
void Reset(uint64_t CurrentSimNanos)¶
This method is used to reset the module.
- Returns
void
-
void UpdateState(uint64_t CurrentSimNanos)¶
update state
- Parameters
currentSimNanos –
-
void linkInStates(DynParamManager &statesIn)¶
class method
-
void registerProperties(DynParamManager &statesIn)¶
register properties
- Parameters
statesIn – simulation states
-
void computeGravityField(Eigen::Vector3d r_cF_N, Eigen::Vector3d rDot_cF_N)¶
Calculate gravitational acceleration of s/c wrt inertial (no central body) or wrt central body
- Parameters
r_cF_N – is position of center of mass of s/c wrt frame it is stored/integrated in in spacecraft
rDot_cF_N – is the derivative of above
-
void updateInertialPosAndVel(Eigen::Vector3d r_BF_N, Eigen::Vector3d rDot_BF_N)¶
Calculate gravitational acceleration of s/c wrt inertial (no central body) or wrt central body
- Parameters
r_BF_N – is position of body frame of s/c wrt frame it is stored/integrated in in spacecraft
rDot_BF_N – is the derivative of above
-
void updateEnergyContributions(Eigen::Vector3d r_CN_N, double &orbPotEnergyContr)¶
— Orbital Potential Energy Contributions
-
void setGravBodies(std::vector<GravBodyData*> gravBodies)¶
class method
-
void addGravBody(GravBodyData *gravBody)¶
class method
-
void prependSpacecraftNameToStates()¶
class method
Public Members
-
std::string vehicleGravityPropName¶
[-] Name of the vehicle mass state
-
std::string systemTimeCorrPropName¶
[-] Name of the correlation between times
-
std::vector<GravBodyData*> gravBodies¶
[-] Vector of bodies we feel gravity from
-
GravBodyData *centralBody¶
Central body.
-
std::string inertialPositionPropName¶
[-] Name of the inertial position property
-
std::string inertialVelocityPropName¶
[-] Name of the inertial velocity property
-
std::string nameOfSpacecraftAttachedTo¶
[-] Name of the s/c this gravity model is attached to
-
BSKLogger bskLogger¶
— BSK Logging
-
Message<SpicePlanetStateMsgPayload> centralBodyOutMsg¶
central planet body state output message
Private Functions
-
Eigen::Vector3d getEulerSteppedGravBodyPosition(GravBodyData *bodyData)¶
class method
compute planet position with Euler integration
- Parameters
bodyData – planet data
-
void writeOutputMessages(uint64_t currentSimNanos)¶
class method
write output message
- Parameters
currentSimNanos –
Private Members
-
Eigen::MatrixXd *gravProperty¶
[-] g_N property for output
-
Eigen::MatrixXd *timeCorr¶
[-] Time correlation property
-
Eigen::MatrixXd *inertialPositionProperty¶
[m] r_N inertial position relative to system spice zeroBase/refBase coordinate frame, property for output.
-
Eigen::MatrixXd *inertialVelocityProperty¶
[m/s] v_N inertial velocity relative to system spice zeroBase/refBase coordinate frame, property for output.
-
GravityEffector()¶