Module: constraintDynamicEffector

Executive Summary

The constraint effector class is used to define a physical connection between two separate spacecraft. This class takes in the connection location on each spacecraft in their respective body frames as well as a relative vector between these two connection points. It also takes in two gain tuning parameters to tailor the stiffness and damping of the connection between the two spacecraft.

The constraint effector class is an instantiation of the dynamic effector abstract class. The included test is validating the interaction between two spacecraft rigid body hubs that are attached through a constraint effector. In this case, two identical spacecraft are connected by a 0.1 meter long arm which is enforced with high stiffness and damping to be virtually rigid.

Detailed Module Description

A constraint effector is a combination of two separate 3-degree-of-freedom holonomic constraints: a direction constraint which enforces the position of the connection point on each spacecraft relative to the other, and an attitude constraint which enforces the attitude of each spacecraft relative to the other. The constraint effector works by correcting violations of these constraints as well as their derivatives.

Mathematical Modeling

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

Note

J. Vaz Carneiro, A. Morell and H. Schaub, “Post-Docking Complex Spacecraft Dynamics Using Baumgarte Stabilization”, AAS/AIAA Astrodynamics Specialist Conference, Big Sky, Montana, Aug. 13-17, 2023

User Guide

This section outlines the steps needed to setup a Constraint Dynamic Effector in Python using Basilisk.

  1. Import the constraintDynamicEffector class:

    from Basilisk.simulation import constraintDynamicEffector
    
  2. Create an instantiation of a holonomic constraint:

    constraintEffector = constraintDynamicEffector.ConstraintDynamicEffector()
    
  3. Define all physical parameters for the constraint:

    constraintEffector.setR_P1B1_B1(r_P1B1_B1)
    constraintEffector.setR_P2B2_B2(r_P2B2_B2)
    constraintEffector.setR_P2P1_B1Init(r_P2P1_B1Init)
    
  4. Define the stiffness and damping of the connection. See the recommended starting gains here:

    constraintEffector.setAlpha(1E2)
    constraintEffector.setBeta(1e3)
    
  5. (Optional) Define exact gains for the direction and attitude constraints separately. These are internally set based on alpha and beta during reset, but can be overridden in this way:

    constraintEffector.setK_d(alpha**2)
    constraintEffector.setC_d(2*beta)
    constraintEffector.setK_a(alpha**2)
    constraintEffector.setC_a(2*beta)
    
  6. Add the effector to both spacecraft:

    scObject1.addDynamicEffector(constraintEffector)
    scObject2.addDynamicEffector(constraintEffector)
    

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

  7. Add the module to the task list:

    Sim.AddModelToTask(TaskName, constraintEffector)
    
  8. Ensure that the dynamics integration is synced between the connected spacecraft. For best results use a variable timestep integrator:

    integratorObject1 = svIntegrators.svIntegratorRKF45(scObject1)
    scObject1.setIntegrator(integratorObject1)
    scObject1.syncDynamicsIntegration(scObject2)
    
  9. (Optional) Retrieve the constraint effector’s physical parameters, gain tuning parameters, or exact gains for the direction and attitude constraints:

    constraintEffector.getR_P1B1_B1(r_P1B1_B1)
    constraintEffector.getR_P2B2_B2(r_P2B2_B2)
    constraintEffector.getR_P2P1_B1Init(r_P2P1_B1Init)
    constraintEffector.getAlpha(1E2)
    constraintEffector.getBeta(1e3)
    constraintEffector.getK_d(alpha**2)
    constraintEffector.getC_d(2*beta)
    constraintEffector.getK_a(alpha**2)
    constraintEffector.getC_a(2*beta)
    
  10. (Optional) Define a input device status message.(1 means constraintEffector is connected.0 means constraintEffector is disconnected). If not set, it defaults to being connected:

    effectorStatusMsgPayload = messaging.DeviceStatusMsgPayload()
    effectorStatusMsgPayload.deviceStatus = 1
    effectorStatusMsg = messaging.DeviceStatusMsg().write(effectorStatusMsgPayload)
    constraintEffector.effectorStatusInMsg.subscribeTo(effectorStatusMsg)
    
  11. (Optional) Setup Low Pass Filtering for the Constraint Forces and Torques acting on the two satellites. Define a positive cut-off frequency wc for the low-pass filter. If not set, defaults to 0:

    constraintEffector.setFilterData(wc)
    
  12. The constraintEffector output message records the raw and filtered constraint forces and torques acting on the two spacecraft using the variable constraintElements.


class ConstraintDynamicEffector : public SysModel, public DynamicEffector
#include <constraintDynamicEffector.h>

constraint dynamic effector class

Public Functions

ConstraintDynamicEffector()

This is the constructor, nothing to report here

~ConstraintDynamicEffector()

This is the destructor, nothing to report here

void Reset(uint64_t CurrentSimNanos)

This method is used to reset the module.

Returns:

void

void linkInStates(DynParamManager &states)

This method allows the constraint effector to have access to the parent states

Parameters:

states – The states to link

Returns:

void

void computeForceTorque(double integTime, double timeStep)

This method computes the forces on torques on each spacecraft body.

Parameters:
  • integTime – Integration time

  • timeStep – Current integration time step used

Returns:

void

void UpdateState(uint64_t CurrentSimNanos)

Update state method

Parameters:

CurrentSimNanos – The current simulation time

Returns:

void

void writeOutputStateMessage(uint64_t CurrentClock)

This method takes the computed constraint force and torque states and outputs them to the m messaging system.

Parameters:

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

Returns:

void

void computeFilteredForce(uint64_t CurrentClock)

Filtering method to calculate filtered Constraint Force

Parameters:

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

Returns:

void

void computeFilteredTorque(uint64_t CurrentClock)

Filtering method to calculate filtered Constraint Torque

Parameters:

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

Returns:

void

void readInputMessage()

This method allows the user to set the status of the constraint dynamic effector

Returns:

void

void setR_P2P1_B1Init(Eigen::Vector3d r_P2P1_B1Init)

setter for r_P2P1_B1Init initial spacecraft separation

void setR_P1B1_B1(Eigen::Vector3d r_P1B1_B1)

setter for r_P1B1_B1 connection point position on spacecraft 1

void setR_P2B2_B2(Eigen::Vector3d r_P2B2_B2)

setter for r_P2B2_B2 connection point position on spacecraft 2

void setAlpha(double alpha)

setter for alpha gain tuning parameter

void setBeta(double beta)

setter for beta gain tuning parameter

void setK_d(double k_d)

setter for k_d gain

void setC_d(double c_d)

setter for c_d gain

void setK_a(double k_a)

setter for k_a gain

void setC_a(double c_a)

setter for c_a gain

void setFilter_Data(double wc, double h, double k)

setter for a,b,s,c,d,e coefficients of low pass filter

This method allows the user to set the cut-off frequency of the low pass filter which is then used to calculate the coefficients for numerical low pass filtering based on a second-order low pass filter design.

Parameters:
  • wc – The cut-off frequency of the low pass filter.

  • h – The constant digital time step.

  • k – The damping coefficient

Returns:

void

inline Eigen::Vector3d getR_P2P1_B1Init() const

getter for r_P2P1_B1Init initial spacecraft separation

inline Eigen::Vector3d getR_P1B1_B1() const

getter for r_P1B1_B1 connection point position on spacecraft 1

inline Eigen::Vector3d getR_P2B2_B2() const

getter for r_P2B2_B2 connection point position on spacecraft 2

inline double getAlpha() const

getter for alpha gain tuning parameter

inline double getBeta() const

getter for beta gain tuning parameter

inline double getK_d() const

getter for k_d gain

inline double getC_d() const

getter for c_d gain

inline double getK_a() const

getter for k_a gain

inline double getC_a() const

getter for c_a gain

Public Members

Message<ConstDynEffectorMsgPayload> constraintElements

output message with constraint force and torque on connected s/c

ReadFunctor<DeviceStatusMsgPayload> effectorStatusInMsg

input message to record device status

uint64_t effectorStatus = 1

internal variable to toggle effector status

Private Members

int scInitCounter = 0

counter to kill simulation if more than two spacecraft initialized

int scID = 1

0,1 alternating spacecraft tracker to output appropriate force/torque

Eigen::Vector3d r_P1B1_B1 = Eigen::Vector3d::Zero()

[m] position vector from spacecraft 1 hub to its connection point P1

Eigen::Vector3d r_P2B2_B2 = Eigen::Vector3d::Zero()

[m] position vector from spacecraft 2 hub to its connection point P2

Eigen::Vector3d r_P2P1_B1Init = Eigen::Vector3d::Zero()

[m] precribed position vector from spacecraft 1 connection point to spacecraft 2 connection point

double alpha = 0.0

Baumgarte stabilization gain tuning variable.

double beta = 0.0

Baumgarte stabilization gain tuning variable.

double k_d = 0.0

direction constraint proportional gain

double c_d = 0.0

direction constraint derivative gain

double k_a = 0.0

attitude constraint proportional gain

double c_a = 0.0

attitude constraint derivative gain

double a = 0.0

coefficient in numerical low pass filter

double b = 0.0

coefficient in numerical low pass filter

double c = 0.0

coefficient in numerical low pass filter

double d = 0.0

coefficient in numerical low pass filter

double e = 0.0

coefficient in numerical low pass filter

double F_mag_tminus2 = 0.0

Magnitude of unfiltered constraint force at t-2 time step.

double F_mag_tminus1 = 0.0

Magnitude of unfiltered constraint force at t-1 time step.

double F_mag_t = 0.0

Magnitude of unfiltered constraint force at t time step.

double F_filtered_mag_t = 0.0

Magnitude of filtered constraint force at t time step.

double F_filtered_mag_tminus1 = 0.0

Magnitude of filtered constraint force at t-1 time step.

double F_filtered_mag_tminus2 = 0.0

Magnitude of filtered constraint force at t-2 time step.

double T1_mag_tminus2 = 0.0

Magnitude of unfiltered constraint torque on s/c 1 at t-2 time step.

double T1_mag_tminus1 = 0.0

Magnitude of unfiltered constraint torque on s/c 1 at t-1 time step.

double T1_mag_t = 0.0

Magnitude of unfiltered constraint torque on s/c 1 at t time step.

double T1_filtered_mag_t = 0.0

Magnitude of filtered constraint torque on s/c 1 at t time step.

double T1_filtered_mag_tminus1 = 0.0

Magnitude of filtered constraint torque on s/c 1 at t-1 time step.

double T1_filtered_mag_tminus2 = 0.0

Magnitude of filtered constraint torque on s/c 1 at t-2 time step.

double T2_mag_tminus2 = 0.0

Magnitude of unfiltered constraint torque on s/c 2 at t-2 time step.

double T2_mag_tminus1 = 0.0

Magnitude of unfiltered constraint torque on s/c 2at t-1 time step.

double T2_mag_t = 0.0

Magnitude of unfiltered constraint torque on s/c 2 at t time step.

double T2_filtered_mag_t = 0.0

Magnitude of filtered constraint torque on s/c 2 at t time step.

double T2_filtered_mag_tminus1 = 0.0

Magnitude of filtered constraint torque on s/c 2 at t-1 time step.

double T2_filtered_mag_tminus2 = 0.0

Magnitude of filtered constraint torque on s/c 2 at t-2 time step.

std::vector<StateData*> hubPosition

[m] parent inertial position vector

std::vector<StateData*> hubVelocity

[m/s] parent inertial velocity vector

std::vector<StateData*> hubSigma

parent attitude Modified Rodrigues Parameters (MRPs)

std::vector<StateData*> hubOmega

[rad/s] parent inertial angular velocity vector

Eigen::Vector3d psi_N = Eigen::Vector3d::Zero()

[m] direction constraint violation in inertial frame

Eigen::Vector3d psiPrime_N = Eigen::Vector3d::Zero()

[m/s] direction rate constraint violation in inertial frame

Eigen::MRPd sigma_B2B1 = Eigen::MRPd::Identity()

attitude constraint violation

Eigen::Vector3d omega_B2B1_B2 = Eigen::Vector3d::Zero()

[rad/s] angular velocity constraint violation in spacecraft 2 body frame

Eigen::Vector3d Fc_N = Eigen::Vector3d::Zero()

[N] force applied on each spacecraft COM in the inertial frame

Eigen::Vector3d T_B2 = Eigen::Vector3d::Zero()

[N-m] torque applied on spacecraft 2 in its body frame

Eigen::Vector3d T_B1 = Eigen::Vector3d::Zero()

[N-m] torque applied on spacecraft 1 in its body frame