Module: spinningBodyTwoDOFStateEffector

Executive Summary

The two-degree-of-freedom spinning body class is an instantiation of the state effector abstract class. The integrated test is validating the interaction between the 2-DoF spinning body module and the rigid body hub that it is attached to. In this case, a 2-DoF spinning body system has two masses and two inertia tensors. The lower axis is attached to the hub and the upper axis is attached to the lower body. This module can represent multiple different effectors, such as a dual-hinged solar array, a control moment gyroscope or a dual-gimbal. A spring and damper can be included in each axis, and an optional motor torque can be applied on each spinning axis.

Nominally, two rigid bodies can be defined, which would represent a chain of 1D rotary joints. However, by setting the mass and the inertia of the lower spinning body to 0, the module can simulate a single rigid component rotating about two axis instead.

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.

Module I/O Messages

Msg Variable Name

Msg Type

Description

spinningBodyOutMsgs

HingedRigidBodyMsgPayload

Output vector of messages containing the spinning body state angle and angle rate

motorTorqueInMsg

ArrayMotorTorqueMsgPayload

(Optional) Input message of the motor torque value

motorLockInMsg

ArrayEffectorLockMsgPayload

(Optional) Input message for locking the axis

spinningBodyRefInMsgs

HingedRigidBodyMsgPayload

(Optional) Input array of messages for prescribing the angles and angle rates

spinningBodyConfigLogOutMsgs

SCStatesMsgPayload

Output vector of messages containing the spinning body inertial position and attitude states

Detailed Module Description

A 2 DoF spinning body has 4 states: theta1, theta2, theta1Dot and theta2Dot.

Mathematical Modeling

See the following conference paper for a detailed description of this model.

Note

J. Vaz Carneiro, C. Allard and H. Schaub, “Rotating Rigid Body Dynamics Architecture for Spacecraft Simulation Software Implementation”, AAS Rocky Mountain GN&C Conference, Breckenridge, CO, Feb. 2–8, 2023

User Guide

This section is to outline the steps needed to setup a Spinning Body 2 DoF State Effector in Python using Basilisk.

  1. Import the spinningBody2DOFStateEffector class:

    from Basilisk.simulation import spinningBody2DOFStateEffector
    
  2. Create an instantiation of a Spinning body:

    spinningBody = spinningBody2DOFStateEffector.SpinningBody2DOFStateEffector()
    
  3. Define all physical parameters for both spinning bodies. For example:

    spinningBody.mass1 = 100.0
    spinningBody.mass2 = 50.0
    spinningBody.IS1PntSc1_S1 = [[100.0, 0.0, 0.0], [0.0, 50.0, 0.0], [0.0, 0.0, 50.0]]
    spinningBody.IS2PntSc2_S2 = [[50.0, 0.0, 0.0], [0.0, 30.0, 0.0], [0.0, 0.0, 40.0]]
    spinningBody.dcm_S10B = [[-1.0, 0.0, 0.0], [0.0, -1.0, 0.0], [0.0, 0.0, 1.0]]
    spinningBody.dcm_S20S1 = [[0.0, -1.0, 0.0], [0.0, .0, -1.0], [1.0, 0.0, 0.0]]
    spinningBody.r_Sc1S1_S1 = [[2.0], [-0.5], [0.0]]
    spinningBody.r_Sc2S2_S2 = [[1.0], [0.0], [-1.0]]
    spinningBody.r_S1B_B = [[-2.0], [0.5], [-1.0]]
    spinningBody.r_S2S1_S1 = [[0.5], [-1.5], [-0.5]]
    spinningBody.s1Hat_S1 = [[0], [0], [1]]
    spinningBody.s2Hat_S2 = [[0], [-1], [0]]
    
  4. (Optional) Define initial conditions of the effector. Default values are zero states:

    spinningBody.theta1Init = 0 * macros.D2R
    spinningBody.theta1DotInit = 0.1 * macros.D2R
    spinningBody.theta2Init = 5 * macros.D2R
    spinningBody.theta2DotInit = -0.5 * macros.D2R
    
  5. (Optional) Define spring and damper coefficients. Default values are zero states:

    spinningBody.k1 = 1.0
    spinningBody.c1 = 0.1
    spinningBody.k2 = 2.0
    spinningBody.c2 = 0.5
    
  6. (Optional) Define a unique name for each state. If you have multiple spinning bodies, they each must have a unique name. If these names are not specified, then the default names are used which are incremented by the effector number:

    spinningBody.nameOfTheta1State = "spinningBodyTheta1"
    spinningBody.nameOfTheta1DotState = "spinningBodyTheta1Dot"
    spinningBody.nameOfTheta2State = "spinningBodyTheta2"
    spinningBody.nameOfTheta2DotState = "spinningBodyTheta2Dot"
    
  7. (Optional) Connect a command torque message:

    cmdArray = messaging.ArrayMotorTorqueMsgPayload()
    cmdArray.motorTorque = [cmdTorque1, cmdTorque2]  # [Nm]
    cmdMsg = messaging.ArrayMotorTorqueMsg().write(cmdArray)
    spinningBody.motorTorqueInMsg.subscribeTo(cmdMsg)
    
  8. (Optional) Connect an axis-locking message (0 means the axis is free to rotate and 1 locks the axis):

    lockArray = messaging.ArrayEffectorLockMsgPayload()
    lockArray.motorTorque = [1, 0]
    lockMsg = messaging.ArrayEffectorLockMsg().write(lockArray)
    spinningBody.motorLockInMsg.subscribeTo(lockMsg)
    
  9. (Optional) Connect angle and angle rate reference messages:

    angle1Ref = messaging.HingedRigidBodyMsgPayload()
    angle1Ref.theta = theta1Ref
    angle1Ref.thetaDot = theta1DotRef
    angle1RefMsg = messaging.HingedRigidBodyMsg().write(angle1Ref)
    spinningBody.spinningBodyRefInMsgs[0].subscribeTo(angle1RefMsg)
    
    angle2Ref = messaging.HingedRigidBodyMsgPayload()
    angle2Ref.theta = theta2Ref
    angle2Ref.thetaDot = theta2DotRef
    angle2RefMsg = messaging.HingedRigidBodyMsg().write(angle2Ref)
    spinningBody.spinningBodyRefInMsgs[1].subscribeTo(angle2RefMsg)
    
  10. The angular states of the body are created using an output vector of messages spinningBodyOutMsgs.

  11. The spinning body config log state output messages is spinningBodyConfigLogOutMsgs.

  12. Add the effector to your spacecraft:

    scObject.addStateEffector(spinningBody)
    

    See Module: spacecraft documentation on how to set up a spacecraft object.

  13. Add the module to the task list:

    unitTestSim.AddModelToTask(unitTaskName, spinningBody)
    

class SpinningBodyTwoDOFStateEffector : public StateEffector, public SysModel
#include <spinningBodyTwoDOFStateEffector.h>

spinning body state effector class

Public Functions

SpinningBodyTwoDOFStateEffector()

&#8212; Contructor

This is the constructor, setting variables to default values

~SpinningBodyTwoDOFStateEffector()

&#8212; Destructor

This is the destructor, nothing to report here

void Reset(uint64_t CurrentClock)

&#8212; Method for reset

This method is used to reset the module.

void writeOutputStateMessages(uint64_t CurrentClock)

&#8212; Method for writing the output messages

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

void UpdateState(uint64_t CurrentSimNanos)

&#8212; Method for updating information

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

void registerStates(DynParamManager &statesIn)

&#8212; Method for registering the SB states

This method allows the SB state effector to register its states: theta and thetaDot with the dynamic parameter manager

void linkInStates(DynParamManager &states)

&#8212; Method for getting access to other states

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

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

&#8212; Method for back-substitution contributions

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

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

&#8212; Method for SB to compute its derivatives

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

void updateEffectorMassProps(double integTime)

&#8212; Method for giving the s/c the HRB mass props and prop rates

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

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

&#8212; Method for computing energy and momentum for SBs

This method is for calculating the contributions of the SB state effector to the energy and momentum of the spacecraft

void prependSpacecraftNameToStates()

Method used for multiple spacecraft.

This method prepends the name of the spacecraft for multi-spacecraft simulations.

void computeSpinningBodyInertialStates()

Method for computing the SB’s states.

This method computes the spinning body states relative to the inertial frame

Public Members

double mass1 = 0.0

[kg] mass of lower spinning body (can be 0)

double mass2 = 1.0

[kg] mass of upper spinning body

double k1 = 0.0

[N-m/rad] torsional spring constant for first rotation axis

double k2 = 0.0

[N-m/rad] torsional spring constant for second rotation axis

double c1 = 0.0

[N-m-s/rad] rotational damping coefficient for first rotation axis

double c2 = 0.0

[N-m-s/rad] rotational damping coefficient for second rotation axis

double theta1Init = 0.0

[rad] initial first axis angle

double theta1DotInit = 0.0

[rad/s] initial first axis angle rate

double theta2Init = 0.0

[rad] initial second axis angle

double theta2DotInit = 0.0

[rad/s] initial second axis angle rate

std::string nameOfTheta1State

&#8212; identifier for the theta1 state data container

std::string nameOfTheta1DotState

&#8212; identifier for the thetaDot1 state data container

std::string nameOfTheta2State

&#8212; identifier for the theta2 state data container

std::string nameOfTheta2DotState

&#8212; identifier for the thetaDot2 state data container

Eigen::Vector3d r_S1B_B = {0.0, 0.0, 0.0}

[m] vector pointing from body frame B origin to lower spinning frame S1 origin in B frame components

Eigen::Vector3d r_S2S1_S1 = {0.0, 0.0, 0.0}

[m] vector pointing from lower spinning frame S1 origin to upper spinning frame S2 origin in S1 frame components

Eigen::Vector3d r_Sc1S1_S1 = {0.0, 0.0, 0.0}

[m] vector pointing from lower spinning frame S1 origin to point Sc1 (center of mass of the lower spinner) in S1 frame components

Eigen::Vector3d r_Sc2S2_S2 = {0.0, 0.0, 0.0}

[m] vector pointing from upper spinning frame S2 origin to point Sc2 (center of mass of the upper spinner) in S2 frame components

Eigen::Vector3d s1Hat_S1 = {1.0, 0.0, 0.0}

&#8212; first spinning axis in S1 frame components.

Eigen::Vector3d s2Hat_S2 = {1.0, 0.0, 0.0}

&#8212; second spinning axis in S2 frame components.

Eigen::Matrix3d IS1PntSc1_S1

[kg-m^2] Inertia of lower spinning body about point Sc1 in S1 frame components

Eigen::Matrix3d IS2PntSc2_S2

[kg-m^2] Inertia of upper spinning body about point Sc2 in S2 frame components

Eigen::Matrix3d dcm_S10B

&#8212; DCM from the body frame to the S10 frame (S1 frame for theta1=0)

Eigen::Matrix3d dcm_S20S1

&#8212; DCM from the S1 frame to the S20 frame (S2 frame for theta2=0)

std::vector<Message<HingedRigidBodyMsgPayload>*> spinningBodyOutMsgs{new Message<HingedRigidBodyMsgPayload>, new Message<HingedRigidBodyMsgPayload>}

vector of state output messages

std::vector<Message<SCStatesMsgPayload>*> spinningBodyConfigLogOutMsgs{new Message<SCStatesMsgPayload>, new Message<SCStatesMsgPayload>}

vector of spinning body state config log messages

ReadFunctor<ArrayMotorTorqueMsgPayload> motorTorqueInMsg

&#8212; (optional) motor torque input message name

ReadFunctor<ArrayEffectorLockMsgPayload> motorLockInMsg

&#8212; (optional) motor lock input message name

std::vector<ReadFunctor<HingedRigidBodyMsgPayload>> spinningBodyRefInMsgs{ReadFunctor<HingedRigidBodyMsgPayload>(), ReadFunctor<HingedRigidBodyMsgPayload>()}

(optional) vector of spinning body reference input messages

Private Members

double u1 = 0.0

[N-m] optional motor torque for first axis

double u2 = 0.0

[N-m] optional motor torque for second axis

int lockFlag1 = 0

[] flag for locking the first rotation axis

int lockFlag2 = 0

[] flag for locking the second rotation axis

double theta1Ref = 0.0

[rad] spinning body reference angle

double theta1DotRef = 0.0

[rad] spinning body reference angle rate

double theta2Ref = 0.0

[rad] spinning body reference angle

double theta2DotRef = 0.0

[rad] spinning body reference angle rate

double mass = 1.0

[kg] mass of the spinner system

Eigen::Matrix<double, 2, 3> ATheta

&#8212; rDDot_BN term for back substitution

Eigen::Matrix<double, 2, 3> BTheta

&#8212; omegaDot_BN term for back substitution

Eigen::Vector2d CTheta = {0.0, 0.0}

&#8212; scalar term for back substitution

Eigen::Vector3d s1Hat_B = {1.0, 0.0, 0.0}

&#8212; first spinning axis in B frame components

Eigen::Vector3d s2Hat_B = {1.0, 0.0, 0.0}

&#8212; second spinning axis in B frame components

Eigen::Vector3d r_Sc1S1_B = {0.0, 0.0, 0.0}

[m] vector pointing from lower spinning frame S1 origin to point Sc1 in B frame components

Eigen::Vector3d r_Sc1B_B = {0.0, 0.0, 0.0}

[m] vector pointing from body frame B origin to point Sc1 in B frame components.

Eigen::Vector3d r_Sc2S2_B = {0.0, 0.0, 0.0}

[m] vector pointing from upper spinning frame S2 origin to point Sc2 in B frame components

Eigen::Vector3d r_S2S1_B = {0.0, 0.0, 0.0}

[m] vector pointing from lower spinning frame S1 origin to upper spinning frame S2 origin in B frame components

Eigen::Vector3d r_Sc2S1_B = {0.0, 0.0, 0.0}

[m] vector pointing from lower spinning frame S1 origin to point Sc2 in B frame components

Eigen::Vector3d r_Sc2B_B = {0.0, 0.0, 0.0}

[m] vector pointing from body frame B origin to point Sc2 in B frame components.

Eigen::Vector3d r_ScB_B = {0.0, 0.0, 0.0}

[m] vector pointing from body frame B origin to point Sc (center of mass of the spinner system) in B frame components.

Eigen::Vector3d rPrime_Sc1S1_B = {0.0, 0.0, 0.0}

[m/s] body frame time derivative of r_Sc1S1_B

Eigen::Vector3d rPrime_Sc2S2_B = {0.0, 0.0, 0.0}

[m/s] body frame time derivative of r_Sc2S2_B

Eigen::Vector3d rPrime_S2S1_B = {0.0, 0.0, 0.0}

[m/s] body frame time derivative of r_S2S1_B

Eigen::Vector3d rPrime_Sc1B_B = {0.0, 0.0, 0.0}

[m/s] body frame time derivative of r_Sc1B_B

Eigen::Vector3d rPrime_Sc2B_B = {0.0, 0.0, 0.0}

[m/s] body frame time derivative of r_Sc2B_B

Eigen::Vector3d rPrime_Sc2S1_B = {0.0, 0.0, 0.0}

[m/s] body frame time derivative of r_Sc2S1_B

Eigen::Vector3d rPrime_ScB_B = {0.0, 0.0, 0.0}

[m/s] body frame time derivative of r_ScB_B

Eigen::Vector3d rDot_Sc1B_B = {0.0, 0.0, 0.0}

[m/s] inertial frame time derivative of r_Sc1B_B

Eigen::Vector3d rDot_Sc2B_B = {0.0, 0.0, 0.0}

[m/s] inertial frame time derivative of r_Sc2B_B

Eigen::Vector3d omega_S1B_B = {0.0, 0.0, 0.0}

[rad/s] angular velocity of the S1 frame wrt the B frame in B frame components

Eigen::Vector3d omega_S2S1_B = {0.0, 0.0, 0.0}

[rad/s] angular velocity of the S2 frame wrt the S1 frame in B frame components

Eigen::Vector3d omega_S2B_B = {0.0, 0.0, 0.0}

[rad/s] angular velocity of the S2 frame wrt the B frame in B frame components

Eigen::Vector3d omega_BN_B = {0.0, 0.0, 0.0}

[rad/s] angular velocity of the B frame wrt the N frame in B frame components

Eigen::Vector3d omega_S1N_B = {0.0, 0.0, 0.0}

[rad/s] angular velocity of the S1 frame wrt the N frame in B frame components

Eigen::Vector3d omega_S2N_B = {0.0, 0.0, 0.0}

[rad/s] angular velocity of the S2 frame wrt the N frame in B frame components

Eigen::MRPd sigma_BN = {0.0, 0.0, 0.0}

&#8212; body frame attitude wrt to the N frame in MRPs

Eigen::Matrix3d rTilde_Sc1B_B

[m] tilde matrix of r_Sc1B_B

Eigen::Matrix3d rTilde_Sc2B_B

[m] tilde matrix of r_Sc2B_B

Eigen::Matrix3d omegaTilde_S1B_B

[rad/s] tilde matrix of omega_S1B_B

Eigen::Matrix3d omegaTilde_S2B_B

[rad/s] tilde matrix of omega_S2B_B

Eigen::Matrix3d omegaTilde_BN_B

[rad/s] tilde matrix of omega_BN_B

Eigen::Matrix3d dcm_BS1

&#8212; DCM from lower spinner frame to body frame

Eigen::Matrix3d dcm_BS2

&#8212; DCM from upper spinner frame to body frame

Eigen::Matrix3d dcm_BN

&#8212; DCM from inertial frame to body frame

Eigen::Matrix3d IS1PntSc1_B

[kg-m^2] inertia of lower spinning body about point Sc1 in B frame components

Eigen::Matrix3d IS2PntSc2_B

[kg-m^2] inertia of upper spinning body about point Sc2 in B frame components

Eigen::Matrix3d IPrimeS1PntSc1_B

[kg-m^2] body frame time derivative of inertia of inertia of lower spinning body about point Sc1 in B frame components

Eigen::Matrix3d IPrimeS2PntSc2_B

[kg-m^2] body frame time derivative of inertia of inertia of upper spinning body about point Sc2 in B frame components

Eigen::Vector3d r_Sc1N_N = {0.0, 0.0, 0.0}

[m] position vector of lower spinning body center of mass Sc1 relative to the inertial frame origin N

Eigen::Vector3d r_Sc2N_N = {0.0, 0.0, 0.0}

[m] position vector of upper spinning body center of mass Sc2 relative to the inertial frame origin N

Eigen::Vector3d v_Sc1N_N = {0.0, 0.0, 0.0}

[m/s] inertial velocity vector of Sc1 relative to inertial frame

Eigen::Vector3d v_Sc2N_N = {0.0, 0.0, 0.0}

[m/s] inertial velocity vector of Sc2 relative to inertial frame

Eigen::Vector3d sigma_S1N = {0.0, 0.0, 0.0}

&#8212; MRP attitude of frame S1 relative to inertial frame

Eigen::Vector3d sigma_S2N = {0.0, 0.0, 0.0}

&#8212; MRP attitude of frame S2 relative to inertial frame

Eigen::Vector3d omega_S1N_S1 = {0.0, 0.0, 0.0}

[rad/s] inertial lower spinning body frame angular velocity vector

Eigen::Vector3d omega_S2N_S2 = {0.0, 0.0, 0.0}

[rad/s] inertial upper spinning body frame angular velocity vector

double theta1 = 0.0

[rad] first axis angle

double theta1Dot = 0.0

[rad/s] first axis angle rate

double theta2 = 0.0

[rad] second axis angle

double theta2Dot = 0.0

[rad/s] second axis angle rate

Eigen::MatrixXd *inertialPositionProperty = nullptr

[m] r_N inertial position relative to system spice zeroBase/refBase

Eigen::MatrixXd *inertialVelocityProperty = nullptr

[m] v_N inertial velocity relative to system spice zeroBase/refBase

StateData *theta1State = nullptr

&#8212; state manager of theta1 for spinning body

StateData *theta1DotState = nullptr

&#8212; state manager of theta1Dot for spinning body

StateData *theta2State = nullptr

&#8212; state manager of theta2 for spinning body

StateData *theta2DotState = nullptr

&#8212; state manager of theta2Dot for spinning body

Private Static Attributes

static uint64_t effectorID = 1

[] ID number of this panel