algebraic_riccati_continuous.h
Go to the documentation of this file.
1 /*********************************************************************
2  *
3  * Software License Agreement
4  *
5  * Copyright (c) 2020,
6  * TU Dortmund - Institute of Control Theory and Systems Engineering.
7  * All rights reserved.
8  *
9  * This program is free software: you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation, either version 3 of the License, or
12  * (at your option) any later version.
13  *
14  * This program is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the GNU General Public License
20  * along with this program. If not, see <https://www.gnu.org/licenses/>.
21  *
22  * Authors: Christoph Rösmann
23  *********************************************************************/
24 
25 #ifndef SRC_NUMERICS_INCLUDE_CORBO_NUMERICS_ALGEBRAIC_RICCATI_CONTINUOUS_H_
26 #define SRC_NUMERICS_INCLUDE_CORBO_NUMERICS_ALGEBRAIC_RICCATI_CONTINUOUS_H_
27 
28 #include <corbo-core/types.h>
29 
30 #include <Eigen/Eigenvalues>
31 
32 namespace corbo {
33 
58 class AlgebraicRiccatiContinuous
59 {
60  public:
84  const Eigen::Ref<const Eigen::MatrixXd>& Q, const Eigen::Ref<const Eigen::MatrixXd>& R, Eigen::MatrixXd& X,
85  Eigen::MatrixXd* G = nullptr);
86 
102 
119 
120  protected:
145  static bool solveRiccatiHamiltonianSchur(const Eigen::MatrixXd& H, Eigen::MatrixXd& X);
146 
147  // matrix sign method, currently not in use and needs to be updated for future use
148  // [1] R. Byers. "Solving the algebraic Riccati equation with the matrix sign function".
149  // Linear ALgebra and its Applications. Volume 85, 267-278, 1987.
150  // static void sovleRiccatiHamiltonianMatrixSign(Eigen::MatrixXd& Z, Eigen::MatrixXd& X);
151 
153  static bool hasNegativeRealPart(const Eigen::Ref<const Eigen::MatrixXd>& block) { return block.eigenvalues()[0].real() < 0; }
154 
157 };
158 
159 } // namespace corbo
160 
161 #endif // SRC_NUMERICS_INCLUDE_CORBO_NUMERICS_ALGEBRAIC_RICCATI_CONTINUOUS_H_
corbo
Definition: communication/include/corbo-communication/utilities.h:37
corbo::AlgebraicRiccatiContinuous::hasRealPartsCloseToZero
static bool hasRealPartsCloseToZero(const Eigen::Ref< const Eigen::MatrixXd > &matrix)
Check if the real parts of the provided matrix are close to zero.
Definition: algebraic_riccati_continuous.cpp:125
corbo::AlgebraicRiccatiContinuous::isNumericallyStable
static bool isNumericallyStable(const Eigen::Ref< const Eigen::MatrixXd > &A, const Eigen::Ref< const Eigen::MatrixXd > &B, const Eigen::Ref< const Eigen::MatrixXd > &Q, const Eigen::Ref< const Eigen::MatrixXd > &R)
Determine if solving the riccati equation seems to be numerically stable.
Definition: algebraic_riccati_continuous.cpp:218
X
#define X
Definition: icosphere.cpp:20
corbo::AlgebraicRiccatiContinuous::isClosedLoopStable
static bool isClosedLoopStable(const Eigen::Ref< const Eigen::MatrixXd > &A, const Eigen::Ref< const Eigen::MatrixXd > &B, const Eigen::Ref< const Eigen::MatrixXd > &G)
Determine if the closed-loop system is stable.
Definition: algebraic_riccati_continuous.cpp:202
matrix
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:102
corbo::AlgebraicRiccatiContinuous::hasNegativeRealPart
static bool hasNegativeRealPart(const Eigen::Ref< const Eigen::MatrixXd > &block)
Predicate for Schur block ordering (see solveRiccatiHamiltonianSchur())
Definition: algebraic_riccati_continuous.h:197
B
MatrixType B(b, *n, *nrhs, *ldb)
block
EIGEN_DOC_BLOCK_ADDONS_NOT_INNER_PANEL EIGEN_DEVICE_FUNC BlockXpr block(Index startRow, Index startCol, Index blockRows, Index blockCols)
This is the const version of block(Index,Index,Index,Index). *‍/.
Definition: BlockMethods.h:64
Eigen::Ref
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:192
types.h
corbo::AlgebraicRiccatiContinuous::solveRiccatiHamiltonianSchur
static bool solveRiccatiHamiltonianSchur(const Eigen::MatrixXd &H, Eigen::MatrixXd &X)
Solve Hamiltonian via Schur method.
Definition: algebraic_riccati_continuous.cpp:94
corbo::AlgebraicRiccatiContinuous::solve
static bool solve(const Eigen::Ref< const Eigen::MatrixXd > &A, const Eigen::Ref< const Eigen::MatrixXd > &B, const Eigen::Ref< const Eigen::MatrixXd > &Q, const Eigen::Ref< const Eigen::MatrixXd > &R, Eigen::MatrixXd &X, Eigen::MatrixXd *G=nullptr)
Solve continuous-time algebraic Riccati equation.
Definition: algebraic_riccati_continuous.cpp:62
A
MatrixType A(a, *n, *n, *lda)


control_box_rst
Author(s): Christoph Rösmann
autogenerated on Wed Mar 2 2022 00:05:36