Static Public Member Functions | List of all members
corbo::LyapunovContinuous Class Reference

Methods for dealing with continuous-time Lyapunov equations. More...

#include <lyapunov_continuous.h>

Static Public Member Functions

static bool hasUniqueSolution (const Eigen::Ref< const Eigen::MatrixXd > &A)
 Determine if the Lyapunov equation exhibits a unique solution. More...
 
static bool solve (const Eigen::Ref< const Eigen::MatrixXd > &A, const Eigen::Ref< const Eigen::MatrixXd > &Q, Eigen::MatrixXd &X)
 Solve continuous-time Lyapunov equation. More...
 

Detailed Description

Methods for dealing with continuous-time Lyapunov equations.

The continuous-time Lyapunov equation is given by

\[ A X + X A^T + Q = 0 \]

If Q is symmetric, solution X is symmetric as well. A and Q must be square.

Remarks
The Lyapunov equation is a special case of the more general Sylvester equation. We still provide this particular class in order to provide a more dedicated implementation in terms of efficiency.
See also
LyapunovDiscrete SylvesterContinuous SylvesterDiscrete AlgebraicRiccatiContinuous AlgebraicRiccatiDiscrete
Author
Christoph Rösmann (chris.nosp@m.toph.nosp@m..roes.nosp@m.mann.nosp@m.@tu-d.nosp@m.ortm.nosp@m.und.d.nosp@m.e)
Todo:

Add support for the generalized Lyapunov equations $ A X E^T + E X A^T + Q = 0 $.

Allow the user to precompute the Schur decomposition for subsequent calls of solve() with varying $ Q $.

Definition at line 56 of file lyapunov_continuous.h.

Member Function Documentation

◆ hasUniqueSolution()

bool corbo::LyapunovContinuous::hasUniqueSolution ( const Eigen::Ref< const Eigen::MatrixXd > &  A)
static

Determine if the Lyapunov equation exhibits a unique solution.

The solution is unique if for Eigenvalues of A it is: $ \lambda_i \neq -\lambda_j, \forall i,j. $

Parameters
[in]AMatrix (might be non-square)
Returns
true if the solution is unique, false otherwise (also if matrix A is not square).

Definition at line 83 of file lyapunov_continuous.cpp.

◆ solve()

bool corbo::LyapunovContinuous::solve ( const Eigen::Ref< const Eigen::MatrixXd > &  A,
const Eigen::Ref< const Eigen::MatrixXd > &  Q,
Eigen::MatrixXd &  X 
)
static

Solve continuous-time Lyapunov equation.

Solve $ A X + X A^T + Q = 0 $ w.r.t. $ X $. If Q is symmetric, solution X is symmetric as well. A and Q must be square. The solution is unique if for Eigenvalues of A it is: $ \lambda_i \neq -\lambda_j, \forall i,j. $

The solution is obtained via Schur decomposition [1-4]. We solve the transformed equation: $ T Y + Y T^T + F = 0 $ with

\[ T = U^T A U, \\ F = U^T Q U, \\ Y = U^T X U. \]

In this transformed representation the system of equations can be solved by backward substitution.

[1] R.H. Bartels, G. W. Stewart. "Algorithm 432: Solution of the matrix equation AX + XB = C". Comm. ACM. 15 (9): 820–826, 1972. [2] V. Simoncini. "Computational Methods for Linear Matrix Equations". SIAM Review 58 (3): 377-441, 2016. [3] https://people.kth.se/~eliasj/NLA/matrixeqs.pdf [4] https://github.com/ajt60gaibb/freeLYAP

Warning
Input matrix sizes are checked only in debug mode and occuring issues immediately break execution
Todo:
Reimplement using Real Schur form in order to improve efficiency (if suitable)
Parameters
[in]ASquare matrix
[in]QSquare matrix with same size as A
[out]XSolution with same size as A and Q (size(X) must not be preallocated a-priori)
Returns
true if a finite solution is found, false otherwise.

Definition at line 37 of file lyapunov_continuous.cpp.


The documentation for this class was generated from the following files:


control_box_rst
Author(s): Christoph Rösmann
autogenerated on Mon Feb 28 2022 22:08:02