ADCS.controller.mtq_w_rw_QP module

class ADCS.controller.mtq_w_rw_QP.MTQ_w_RW_QP(est_sat, p_gain, d_gain, c_gain, h_target=array([0., 0., 0.]))[source]

Bases: MTQ_w_RW_LP

Quadratic–Programming–Based Torque Allocation for Mixed RW–MTQ ADCS.

This controller allocates actuator commands for a mixed-actuation attitude control system composed of reaction wheels (RWs) and magnetorquers (MTQs). It computes bounded actuator inputs that best reproduce a desired body torque in a least-squares sense, while strictly respecting actuator saturation limits.

The class extends MTQ_w_RW_LP and replaces the LP-style directional prioritization with a bounded least-squares formulation solved via scipy.optimize.lsq_linear().

Actuator Model

Let the desired body torque be \(\boldsymbol{\tau}_{\mathrm{des}} \in \mathbb{R}^3\). Stack RW and MTQ commands into a single vector

\[\begin{split}\boldsymbol{u} = \begin{bmatrix} \boldsymbol{u}_{\mathrm{rw}} \\ \boldsymbol{u}_{\mathrm{mtq}} \end{bmatrix} \in \mathbb{R}^{n_{\mathrm{rw}} + n_{\mathrm{mtq}}}.\end{split}\]

Reaction wheels generate torque along their spin axes. With RW unit axes stacked as columns in \(A_{\mathrm{rw}} \in \mathbb{R}^{3 \times n_{\mathrm{rw}}}\), the RW torque is

\[\boldsymbol{\tau}_{\mathrm{rw}} = A_{\mathrm{rw}} \, \boldsymbol{u}_{\mathrm{rw}}.\]

Magnetorquers generate a magnetic dipole \(\boldsymbol{m}\) that produces torque via the cross product with the geomagnetic field \(\boldsymbol{B}\):

\[\boldsymbol{\tau}_{\mathrm{mtq}} = \boldsymbol{m} \times \boldsymbol{B} = -[\boldsymbol{B}]_{\times}\,\boldsymbol{m}.\]

If MTQ axes are stacked in \(A_{\mathrm{mtq}} \in \mathbb{R}^{3 \times n_{\mathrm{mtq}}}\), and \(\boldsymbol{u}_{\mathrm{mtq}}\) scales dipole moments along those axes, then

\[\boldsymbol{\tau}_{\mathrm{mtq}} = \left(-[\boldsymbol{B}]_{\times} A_{\mathrm{mtq}}\right)\boldsymbol{u}_{\mathrm{mtq}}.\]

The combined mapping becomes

\[A_{\mathrm{tot}} = \begin{bmatrix} A_{\mathrm{rw}} & -[\boldsymbol{B}]_{\times} A_{\mathrm{mtq}} \end{bmatrix}, \qquad \boldsymbol{\tau}_{\mathrm{ach}} = A_{\mathrm{tot}} \boldsymbol{u}.\]

Bounded Least Squares Allocation

The allocation is formulated as a box-constrained least-squares problem:

\[\min_{\boldsymbol{u}} \;\frac{1}{2}\,\left\|A_{\mathrm{tot}} \boldsymbol{u} - \boldsymbol{\tau}_{\mathrm{des}}\right\|_2^2\]

subject to actuator limits

\[-u_{i,\max} \le u_i \le u_{i,\max}, \qquad \forall i.\]

This approach minimizes Euclidean torque error and naturally returns the closest feasible torque when constraints or MTQ underactuation prevent exact tracking (e.g., torque requests parallel to \(\boldsymbol{B}\)).

Post-hoc Effectiveness Metric

After solving for \(\boldsymbol{u}^*\), the achieved torque is projected onto the desired torque direction to compute an effectiveness scalar \(\alpha\):

\[\hat{\boldsymbol{\tau}} = \frac{\boldsymbol{\tau}_{\mathrm{des}}}{\|\boldsymbol{\tau}_{\mathrm{des}}\|}, \qquad \alpha = \frac{\boldsymbol{\tau}_{\mathrm{ach}} \cdot \hat{\boldsymbol{\tau}}}{\|\boldsymbol{\tau}_{\mathrm{des}}\|}.\]

Values near \(1\) indicate the request is feasible; smaller values indicate saturation or geometric infeasibility.

param est_sat:

Estimated satellite model providing actuator instances and configuration, via EstimatedSatellite.

type est_sat:

EstimatedSatellite

param p_gain:

Proportional gain used by the parent controller law.

type p_gain:

float

param d_gain:

Derivative gain used by the parent controller law.

type d_gain:

float

param c_gain:

Additional control gain used by the parent controller law.

type c_gain:

float

param h_target:

Target wheel momentum vector used by the parent mixed-actuation controller.

type h_target:

numpy.ndarray | list

allocate_max_torque_in_direction(tau_des, b_body, est_sat)[source]

Allocate bounded RW and MTQ commands that best achieve a desired body torque.

The method constructs the combined actuator mapping

\[A_{\mathrm{tot}} = \begin{bmatrix} A_{\mathrm{rw}} & -[\boldsymbol{B}]_{\times} A_{\mathrm{mtq}} \end{bmatrix}\]

using RW axes from RW actuators and MTQ axes from MTQ actuators contained in EstimatedSatellite.

It then solves the box-constrained least-squares problem

\[\min_{\boldsymbol{u}} \;\frac{1}{2}\,\left\|A_{\mathrm{tot}} \boldsymbol{u} - \boldsymbol{\tau}_{\mathrm{des}}\right\|_2^2, \qquad -u_{i,\max} \le u_i \le u_{i,\max}\]

where \(u_{i,\max}\) are the per-actuator saturation limits (RW torque bounds and MTQ dipole-moment bounds, as represented in the actuator objects). The solver used is scipy.optimize.lsq_linear(), which provides a robust Trust Region Reflective method suitable for rank-deficient cases common in magnetic control.

For the solved command \(\boldsymbol{u}^*\), the achieved torque is

\[\boldsymbol{\tau}_{\mathrm{ach}} = A_{\mathrm{tot}} \boldsymbol{u}^*.\]

A post-hoc effectiveness scalar \(\alpha\) is computed as the achieved torque component along the desired direction, normalized by requested magnitude:

\[\hat{\boldsymbol{\tau}} = \frac{\boldsymbol{\tau}_{\mathrm{des}}}{\|\boldsymbol{\tau}_{\mathrm{des}}\|}, \qquad \alpha = \max\!\left(0,\frac{\boldsymbol{\tau}_{\mathrm{ach}} \cdot \hat{\boldsymbol{\tau}}}{\|\boldsymbol{\tau}_{\mathrm{des}}\|}\right).\]

Edge Cases

If \(\|\boldsymbol{\tau}_{\mathrm{des}}\|\) is near zero, the method returns zero commands for all actuators and \(\alpha = 1\) (trivially satisfied request). If there are no actuators, it returns empty command vectors and \(\alpha = 0\).

param tau_des:

Desired body-frame torque vector \(\boldsymbol{\tau}_{\mathrm{des}}\) in N·m.

type tau_des:

numpy.ndarray

param b_body:

Body-frame geomagnetic field vector \(\boldsymbol{B}\) in tesla.

type b_body:

numpy.ndarray

param est_sat:

Estimated satellite model providing actuator instances and their limits.

type est_sat:

EstimatedSatellite

return:

Tuple of RW commands, MTQ commands, and effectiveness scalar \(\alpha\). The RW command vector has length equal to the number of RW actuators. The MTQ command vector has length equal to the number of MTQ actuators.

rtype:

tuple[numpy.ndarray, numpy.ndarray, float]

Parameters:
Return type:

tuple[ndarray, ndarray, float]

Parameters:
  • est_sat (EstimatedSatellite)

  • p_gain (float)

  • d_gain (float)

  • c_gain (float)

  • h_target (ndarray | list)