Module: linearTimeInvariantSystem
Executive Summary
The LinearTimeInvariantSystem module provides a reusable, stateful continuous-time linear time-invariant (LTI) model of the form \(\dot{\mathbf{x}} = \mathbf{A}\mathbf{x} + \mathbf{B}\mathbf{u}\) and \(\mathbf{y} = \mathbf{C}\mathbf{x} + \mathbf{D}\mathbf{u}\). It is a base class intended to be subclassed with application-specific message mappings for the input vector \(\mathbf{u}\) and output vector \(\mathbf{y}\).
Module Description
The base model computes
where \(\mathbf{x} \in \mathbb{R}^{n}\), \(\mathbf{u} \in \mathbb{R}^{m}\), and \(\mathbf{y} \in \mathbb{R}^{p}\).
The class registers a dynamic state block named "x" through the StatefulSysModel interface and initializes it to zero when \(n>0\). During each update step, subclasses provide \(\mathbf{u}\) through readInput(CurrentSimNanos), and the base class computes \(\dot{\mathbf{x}}\) and \(\mathbf{y}\) before delegating output publication to writeOutput(CurrentSimNanos, y).
If no state is registered (for example, \(n=0\)), the module evaluates only the direct-feedthrough term \(\mathbf{y} = \mathbf{D}\mathbf{u}\) when dimensions permit.
The module includes convenience helpers for standard second-order models:
configureSecondOrder(wn, zeta, k)for SISO systems with transfer function\[G(s) = \frac{k\omega_n^2}{s^2 + 2\zeta\omega_n s + \omega_n^2}\]configureSecondOrder(wn, zeta, k)(vector overload) for decoupled MIMO second-order channels with block-diagonal dynamics.
Message Interfaces
This base class does not define concrete simulation messages.
Subclasses must provide message mappings by overriding:
readInput(CurrentSimNanos)to construct \(\mathbf{u}\) from input messageswriteOutput(CurrentSimNanos, y)to publish \(\mathbf{y}\) to output messages
Module Assumptions and Limitations
The system is continuous-time and integrated by the simulation framework.
Matrix dimensions must be consistent with inferred or overridden input/state/output sizes; inconsistencies are reported during
Reset.When \(n=0\), only direct feedthrough \(\mathbf{y}=\mathbf{D}\mathbf{u}\) is evaluated.
Size inference is matrix-driven in the base class unless overridden by derived classes.
Testing
Unit tests are in _UnitTest/test_linearTimeInvariantSystem.py and are skipped when Basilisk is compiled without MuJoCo support (--mujoco). Coverage includes first-order and second-order response checks for C++ and Python-subclass usage.
-
class LinearTimeInvariantSystem : public StatefulSysModel
- #include <linearTimeInvariantSystem.h>
Linear time invariant system of the form x_dot = A x + B u, y = C x + D u.
* x_dot = A x + B u * y = C x + D u *
Subclasses must implement the input and output mapping methods.
Subclassed by ForceAtSiteLTI, OrbitalElementControl, SingleActuatorLTI
State access
Accessors for the internal state vector.
-
BSKLogger bskLogger
Logger used to report errors and informational messages.
-
StateData *xState = nullptr
Pointer to the registered state data for the state vector x.
This pointer is set during registerStates and is used to read and update the state and its derivative.
-
Eigen::MatrixXd B
Input matrix B.
Defines the input contribution to the dynamics x_dot = A x + B u.
-
Eigen::MatrixXd D
Feedthrough matrix D.
Maps the input vector directly to the output y = C x + D u.
-
Eigen::VectorXd getX()
Get the current state vector x.
If the state has not been registered yet, an empty vector is returned and an error is logged.
- Returns:
Copy of the current state vector.
-
void setX(const Eigen::VectorXd &xin)
Set the current state vector x.
If the state has been registered and the supplied vector has the same dimension as the stored state, the state is updated. Otherwise an error is logged and the call is ignored.
- Parameters:
xin – New state vector.
-
virtual size_t getInputSize() const
Get the dimension of the input vector u.
The size is inferred from the B or D matrices when they are set.
- Returns:
Size of the input vector.
-
virtual size_t getStateSize() const
Get the dimension of the state vector x.
The size is inferred from the A, B, or C matrices when they are set.
- Returns:
Size of the state vector.
-
virtual size_t getOutputSize() const
Get the dimension of the output vector y.
The size is inferred from the C or D matrices when they are set.
- Returns:
Size of the output vector.
-
virtual Eigen::VectorXd readInput(uint64_t CurrentSimNanos) = 0
Read the current input vector u from the simulation.
Subclasses must implement this method to construct the input vector of size getInputSize() from the relevant input messages at the given simulation time.
- Parameters:
CurrentSimNanos – Current simulation time in nanoseconds.
- Returns:
Input vector u.
-
virtual void writeOutput(uint64_t CurrentSimNanos, const Eigen::VectorXd &y) = 0
Write the current output vector y to the simulation.
Subclasses must implement this method to map the output vector of size getOutputSize() to the appropriate output messages at the given simulation time.
- Parameters:
CurrentSimNanos – Current simulation time in nanoseconds.
y – Output vector to write.
Parameter getters and setters
Accessors for the system matrices A, B, C, and D.
-
const Eigen::MatrixXd &getA() const
Get the state matrix A.
- Returns:
Constant reference to the A matrix.
-
void setA(const Eigen::MatrixXd &Ain)
Set the state matrix A.
This matrix defines the homogeneous term in the state equation
* x_dot = A x + B u *
- Parameters:
Ain – New A matrix.
-
const Eigen::MatrixXd &getB() const
Get the input matrix B.
- Returns:
Constant reference to the B matrix.
-
void setB(const Eigen::MatrixXd &Bin)
Set the input matrix B.
This matrix defines the contribution of the input vector u to the state equation
* x_dot = A x + B u *
- Parameters:
Bin – New B matrix.
-
const Eigen::MatrixXd &getC() const
Get the output matrix C.
- Returns:
Constant reference to the C matrix.
-
void setC(const Eigen::MatrixXd &Cin)
Set the output matrix C.
This matrix maps the state vector x to the output equation
* y = C x + D u *
- Parameters:
Cin – New C matrix.
-
const Eigen::MatrixXd &getD() const
Get the feedthrough matrix D.
- Returns:
Constant reference to the D matrix.
-
void setD(const Eigen::MatrixXd &Din)
Set the feedthrough matrix D.
This matrix maps the input vector u directly to the output equation
* y = C x + D u *
- Parameters:
Din – New D matrix.
-
void configureSecondOrder(double wn, double zeta, double k = 1.0)
Configure the A, B, C, D matrices as a standard SISO second order model.
The resulting transfer function is:
* G(s) = k * wn^2 / (s^2 + 2*zeta*wn*s + wn^2) *
State space realization:
* x = [position, velocity]^T * x_dot = [ 0 1 ] x + [0 ] u * [ -wn^2 -2*zeta*wn ] [k * wn^2] * y = [1 0] x *
- Parameters:
wn – Natural frequency in rad/s. Must be positive.
zeta – Damping ratio (dimensionless). Must be non negative.
k – Static gain. Default is unity.
-
void configureSecondOrder(const Eigen::VectorXd &wn, const Eigen::VectorXd &zeta, const Eigen::VectorXd &k)
Configure A, B, C, D as a decoupled MIMO second order model.
Each element of the input vectors defines one independent SISO channel:
* G_i(s) = k_i * wn_i^2 / (s^2 + 2*zeta_i*wn_i*s + wn_i^2) *
The full state vector is
* x = [pos_0, vel_0, pos_1, vel_1, ..., pos_{n-1}, vel_{n-1}]^T *
For n channels, the resulting dimensions are:
* A: (2n x 2n), B: (2n x n), C: (n x 2n), D: (n x n). *
- Parameters:
wn – Vector of natural frequencies in rad/s. Each wn(i) must be positive.
zeta – Vector of damping ratios (dimensionless). Each zeta(i) must be non negative.
k – Vector of static gains. Must have the same length as wn and zeta.
Public Functions
-
LinearTimeInvariantSystem() = default
Default constructor.
-
void registerStates(DynParamRegisterer registerer) override
Register the state vector with the dynamics framework.
If the state dimension is greater than zero, this function allocates a state block of size getStateSize() and initializes it to zero.
- Parameters:
registerer – State registration helper provided by the framework.
-
void UpdateState(uint64_t CurrentSimNanos) override
Compute the state derivative and system output at the current time.
This method reads the current input vector using readInput, computes the state derivative x_dot = A x + B u and the output y = C x + D u, and writes the output using writeOutput. If no state is registered, only the direct feedthrough term y = D u is evaluated if applicable.
- Parameters:
CurrentSimNanos – Current simulation time in nanoseconds.
-
virtual void Reset(uint64_t CurrentSimNanos) override
Perform consistency checks and one time initialization.
This method validates that the dimensions of the system matrices A, B, C, and D are consistent with getStateSize, getInputSize, and getOutputSize. Any inconsistency is reported through the logger.
- Parameters:
CurrentSimNanos – Current simulation time in nanoseconds.
-
BSKLogger bskLogger