sylvester_continuous.cpp
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 
26 
29 
30 #include <Eigen/Core>
31 #include <Eigen/Dense>
32 #include <Eigen/Eigenvalues>
33 #include <algorithm>
34 
35 namespace corbo {
36 
38  const Eigen::Ref<const Eigen::MatrixXd>& C, Eigen::MatrixXd& X)
39 {
40  assert(is_square(A));
41  assert(is_square(B));
42  assert(A.rows() == C.rows());
43  assert(B.rows() == C.cols());
44 
47 
48  if (A_schur.info() == Eigen::NoConvergence) return false;
49  if (B_schur.info() == Eigen::NoConvergence) return false;
50 
51  Eigen::MatrixXcd T_A = A_schur.matrixT(); // create copy since we need to modify the diagonal later
52  const Eigen::MatrixXcd& U_A = A_schur.matrixU();
53  const Eigen::MatrixXcd& T_B = B_schur.matrixT();
54  const Eigen::MatrixXcd& U_B = B_schur.matrixU();
55 
56  const int n = C.rows();
57  const int m = C.cols();
58 
59  // transform rhs
60  Eigen::MatrixXcd F = U_A.adjoint() * C * U_B;
61 
62  Eigen::MatrixXcd Y(n, m);
63  Y.setZero();
64 
66 
67  // solve T_A Y + Y T_B + F = 0
68  // construct transformed solution by forward substitution
69  const int mm = m - 1;
70  for (int k = 0; k < m; ++k)
71  {
72  Eigen::VectorXcd rhs = F.col(k) + Y * T_B.col(k);
73  // move non-zero diagonal element to lhs
74  const std::complex<double> offset = T_B(k, k);
75  T_A.diagonal().array() += offset;
76  // compute SVD
77  // TODO(roesmann): here is place for improving efficiency
78  // and robustness by avoiding the SVD (see references).
79  T_A_svd.compute(T_A);
80  // solve linear system
81  Y.col(k) = T_A_svd.solve(-rhs);
82  // revert offset
83  if (k < mm) T_A.diagonal().array() -= offset;
84  }
85 
86  // transform back
87  X = (U_A * Y * U_B.adjoint()).real();
88 
89  return true;
90 }
91 
93 {
94  if (!is_square(A) || !is_square(B)) return false;
95 
96  Eigen::VectorXcd eigvals_A = A.eigenvalues();
97  Eigen::VectorXcd eigvals_B = B.eigenvalues();
98 
99  for (int i = 0; i < eigvals_A.size(); ++i)
100  {
101  for (int j = 0; j < eigvals_B.size(); ++j)
102  {
103  if (approx_equal(eigvals_A[i].real(), -eigvals_B[j].real()) && approx_equal(eigvals_A[i].imag(), -eigvals_B[j].imag())) return false;
104  }
105  }
106  return true;
107 }
108 
109 } // namespace corbo
corbo::SylvesterContinuous::solve
static bool solve(const Eigen::Ref< const Eigen::MatrixXd > &A, const Eigen::Ref< const Eigen::MatrixXd > &B, const Eigen::Ref< const Eigen::MatrixXd > &C, Eigen::MatrixXd &X)
Solve continuous-time Sylvester equation.
Definition: sylvester_continuous.cpp:59
corbo
Definition: communication/include/corbo-communication/utilities.h:37
corbo::SylvesterContinuous::hasUniqueSolution
static bool hasUniqueSolution(const Eigen::Ref< const Eigen::MatrixXd > &A, const Eigen::Ref< const Eigen::MatrixXd > &B)
Determine if the Sylvester equation exhibits a unique solution.
Definition: sylvester_continuous.cpp:114
real
float real
Definition: datatypes.h:10
X
#define X
Definition: icosphere.cpp:20
value_comparison.h
sylvester_continuous.h
Eigen::NoConvergence
@ NoConvergence
Definition: Constants.h:436
corbo::is_square
bool is_square(const Eigen::MatrixBase< Derived > &matrix)
Determine if a given matrix is square.
Definition: matrix_utilities.h:65
B
MatrixType B(b, *n, *nrhs, *ldb)
Eigen::ComputeThinU
@ ComputeThinU
Definition: Constants.h:385
matrix_utilities.h
imag
const EIGEN_DEVICE_FUNC ImagReturnType imag() const
Definition: CommonCwiseUnaryOps.h:95
Eigen::JacobiSVD
Two-sided Jacobi SVD decomposition of a rectangular matrix.
Definition: ForwardDeclarations.h:258
Eigen::ComputeThinV
@ ComputeThinV
Definition: Constants.h:389
Eigen::Ref
A matrix or vector expression mapping an existing expression.
Definition: Ref.h:192
n
PlainMatrixType mat * n
Definition: eigenvalues.cpp:41
A
MatrixType A(a, *n, *n, *lda)
corbo::approx_equal
bool approx_equal(double a, double b, double epsilon=1e-6)
Check if two doubles are appoximately equal.
Definition: value_comparison.h:78
Eigen::ComplexSchur
Performs a complex Schur decomposition of a real or complex square matrix.
Definition: ComplexSchur.h:51


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