Module: nHingedRigidBodyStateEffector

Executive Summary

This class is an instantiation of the stateEffector class and is a N-hinged rigid body effector. This effector is a rigid body attached to the hub through a torsional spring and damper that approximates a flexible appendage. See Allard, Schaub, and Piggott paper: General Hinged Solar Panel Dynamics Approximating First-Order Spacecraft Flexing for a detailed description of this model. A hinged rigid body has 2 states: theta and thetaDot

The module PDF Description contains further information on this module’s function, how to run it, as well as testing.

Message Connection Descriptions

This state effector has no input or output messages.


struct HingedPanel
#include <nHingedRigidBodyStateEffector.h>

Struct containing all the panel variables. All members are public by default so they can be changed by methods of the N_hingedRigidBodyStateEffector class.

Public Members

double mass = 1.0

[kg] mass of hinged panel

double d = 1.0

[m] distance from hinge point to hinged rigid body center of mass

double k = 1.0

[N-m/rad] torsional spring constant of hinge

double c = 0.0

[N-m-s/rad] rotational damping coefficient of hinge

double thetaInit = 0.0

[rad] Initial hinged rigid body angle

double thetaDotInit = 0.0

[rad/s] Initial hinged rigid body angle rate

Eigen::Matrix3d IPntS_S

[kg-m^2] Inertia of hinged rigid body about point S in S frame components

double theta = 0.0

[rad] hinged rigid body angle

double theta_0 = 0.0

[rad] hinged rigid body rest angle

double thetaDot = 0.0

[rad/s] hinged rigid body angle rate

Eigen::Matrix3d dcm_SS_prev

&#8212; DCM from previous S frame to current S frame

Eigen::Matrix3d dcm_SB

&#8212; DCM from body to S frame

Eigen::Vector3d omega_BN_S

[rad/s] omega_BN in S frame components

Eigen::Vector3d omega_SB_B

[rad/s] omega_SB in B frame components

Eigen::Vector3d sHat1_B

&#8212; unit direction vector for the first axis of the S frame

Eigen::Vector3d sHat2_B

&#8212; unit direction vector for the second axis of the S frame

Eigen::Vector3d sHat3_B

&#8212; unit direction vector for the third axis of the S frame

Eigen::Vector3d r_SB_B

&#8212; Vector pointing from B to CoM of hinged rigid body in B frame components

Eigen::Matrix3d rTilde_SB_B

&#8212; Tilde matrix of rSB_B

Eigen::Vector3d rPrime_SB_B

[m/s] Body time derivative of rSB_B

Eigen::Matrix3d rPrimeTilde_SB_B

&#8212; Tilde matrix of rPrime_SB_B

Eigen::Matrix3d ISPrimePntS_B

[kg-m^2/s] time body derivative IPntS in body frame components

class NHingedRigidBodyStateEffector : public StateEffector, public SysModel
#include <nHingedRigidBodyStateEffector.h>

NHingedRigidBodyStateEffector class.

Public Functions

inline void addHingedPanel(HingedPanel NewPanel)

class method

NHingedRigidBodyStateEffector()

&#8212; Contructor

This is the constructor, setting variables to default values

~NHingedRigidBodyStateEffector()

&#8212; Destructor

This is the destructor, nothing to report here

double HeaviFunc(double cond)

&#8212; Heaviside function used for matrix contributions

  • Method for defining the Heaviside function for the EOMs */

void WriteOutputMessages(uint64_t CurrentClock)

This method takes the computed theta states and outputs them to the messaging system.

Parameters:

CurrentClock – The current simulation time (used for time stamping)

Returns:

void

void UpdateState(uint64_t CurrentSimNanos)

This method is used so that the simulation will ask HRB to update messages.

Parameters:

CurrentSimNanos – The current simulation time in nanoseconds

Returns:

void

void registerStates(DynParamManager &statesIn)

&#8212; Method for registering the HRB states

This method allows the HRB state effector to register its states: theta and thetaDot with the dyn param manager

void linkInStates(DynParamManager &states)

&#8212; Method for getting access to other states

This method allows the HRB state effector to have access to the hub states and gravity

void updateEffectorMassProps(double integTime)

&#8212; Method for stateEffector to give mass contributions

This method allows the HRB state effector to provide its contributions to the mass props and mass prop rates of the spacecraft

void updateContributions(double integTime, BackSubMatrices &backSubContr, Eigen::Vector3d sigma_BN, Eigen::Vector3d omega_BN_B, Eigen::Vector3d g_N)

&#8212; Back-sub contributions

This method allows the HRB state effector to give its contributions to the matrices needed for the back-sub method

void updateEnergyMomContributions(double integTime, Eigen::Vector3d &rotAngMomPntCContr_B, double &rotEnergyContr, Eigen::Vector3d omega_BN_B)

&#8212; Energy and momentum calculations

This method is for calculating the contributions of the HRB state effector to the energy and momentum of the s/c

void computeDerivatives(double integTime, Eigen::Vector3d rDDot_BN_N, Eigen::Vector3d omegaDot_BN_B, Eigen::Vector3d sigma_BN)

&#8212; Method for each stateEffector to calculate derivatives

This method is used to find the derivatives for the HRB stateEffector: thetaDDot and the kinematic derivative

void readInputMessages()

&#8212; method to read input messages

This method reads necessary input messages

Returns:

void

Public Members

std::string nameOfThetaState

&#8212; Identifier for the theta state data container

std::string nameOfThetaDotState

&#8212; Identifier for the thetaDot state data container

Eigen::Vector3d r_HB_B

[m] vector pointing from body frame origin to the first Hinge location

Eigen::Matrix3d rTilde_HB_B

&#8212; Tilde matrix of rHB_B

Eigen::Matrix3d dcm_HB

&#8212; DCM from body frame to hinge frame

BSKLogger bskLogger

&#8212; BSK Logging

Private Members

double totalMass

[kg] Total mass of effector

StateData *thetaState

&#8212; state manager of theta for hinged rigid body

StateData *thetaDotState

&#8212; state manager of thetaDot for hinged rigid body

std::vector<HingedPanel> PanelVec

&#8212; vector containing all the info on the different panels

Eigen::MatrixXd matrixADHRB

[-] term needed for back substitution

Eigen::MatrixXd matrixEDHRB

[-] term needed for back substitution

Eigen::MatrixXd matrixFDHRB

[-] term needed for back substitution

Eigen::MatrixXd matrixGDHRB

[-] term needed for back substitution

Eigen::MatrixXd matrixHDHRB

[-] term needed for back substitution

Eigen::MatrixXd matrixKDHRB

[-] term needed for back substitution

Eigen::MatrixXd matrixLDHRB

[-] term needed for back substitution

Eigen::MatrixXd matrixMDHRB

[-] term needed for back substitution

Eigen::VectorXd vectorVDHRB

[-] term needed for back substitution

Eigen::Vector3d aTheta

&#8212; term needed for back substitution

Eigen::Vector3d bTheta

&#8212; term needed for back substitution

Eigen::Vector3d omegaLoc_BN_B

[rad/s] local copy of omegaBN

Eigen::Matrix3d omegaTildeLoc_BN_B

&#8212; tilde matrix of omegaBN

Eigen::MatrixXd *g_N

[m/s^2] Gravitational acceleration in N frame components

Private Static Attributes

static uint64_t effectorID = 1

[] ID number of this panel