10#ifndef EIGEN_MATRIX_SQUARE_ROOT
11#define EIGEN_MATRIX_SQUARE_ROOT
14#include "./InternalHeaderCheck.h"
22template <
typename MatrixType,
typename ResultType>
23void matrix_sqrt_quasi_triangular_2x2_diagonal_block(
const MatrixType& T,
Index i, ResultType& sqrtT) {
26 typedef typename traits<MatrixType>::Scalar Scalar;
27 Matrix<Scalar, 2, 2> block = T.template block<2, 2>(i, i);
28 EigenSolver<Matrix<Scalar, 2, 2> > es(block);
29 sqrtT.template block<2, 2>(i, i) =
30 (es.eigenvectors() * es.eigenvalues().cwiseSqrt().asDiagonal() * es.eigenvectors().inverse()).real();
36template <
typename MatrixType,
typename ResultType>
37void matrix_sqrt_quasi_triangular_1x1_off_diagonal_block(
const MatrixType& T,
Index i,
Index j, ResultType& sqrtT) {
38 typedef typename traits<MatrixType>::Scalar Scalar;
39 Scalar tmp = (sqrtT.row(i).segment(i + 1, j - i - 1) * sqrtT.col(j).segment(i + 1, j - i - 1)).value();
40 sqrtT.coeffRef(i, j) = (T.coeff(i, j) - tmp) / (sqrtT.coeff(i, i) + sqrtT.coeff(j, j));
44template <
typename MatrixType,
typename ResultType>
45void matrix_sqrt_quasi_triangular_1x2_off_diagonal_block(
const MatrixType& T,
Index i,
Index j, ResultType& sqrtT) {
46 typedef typename traits<MatrixType>::Scalar Scalar;
47 Matrix<Scalar, 1, 2> rhs = T.template block<1, 2>(i, j);
48 if (j - i > 1) rhs -= sqrtT.block(i, i + 1, 1, j - i - 1) * sqrtT.block(i + 1, j, j - i - 1, 2);
49 Matrix<Scalar, 2, 2> A = sqrtT.coeff(i, i) * Matrix<Scalar, 2, 2>::Identity();
50 A += sqrtT.template block<2, 2>(j, j).transpose();
51 sqrtT.template block<1, 2>(i, j).transpose() = A.fullPivLu().solve(rhs.transpose());
55template <
typename MatrixType,
typename ResultType>
56void matrix_sqrt_quasi_triangular_2x1_off_diagonal_block(
const MatrixType& T,
Index i,
Index j, ResultType& sqrtT) {
57 typedef typename traits<MatrixType>::Scalar Scalar;
58 Matrix<Scalar, 2, 1> rhs = T.template block<2, 1>(i, j);
59 if (j - i > 2) rhs -= sqrtT.block(i, i + 2, 2, j - i - 2) * sqrtT.block(i + 2, j, j - i - 2, 1);
60 Matrix<Scalar, 2, 2> A = sqrtT.coeff(j, j) * Matrix<Scalar, 2, 2>::Identity();
61 A += sqrtT.template block<2, 2>(i, i);
62 sqrtT.template block<2, 1>(i, j) = A.fullPivLu().solve(rhs);
66template <
typename MatrixType>
67void matrix_sqrt_quasi_triangular_solve_auxiliary_equation(MatrixType& X,
const MatrixType& A,
const MatrixType& B,
68 const MatrixType& C) {
69 typedef typename traits<MatrixType>::Scalar Scalar;
70 Matrix<Scalar, 4, 4> coeffMatrix = Matrix<Scalar, 4, 4>::Zero();
71 coeffMatrix.coeffRef(0, 0) = A.coeff(0, 0) + B.coeff(0, 0);
72 coeffMatrix.coeffRef(1, 1) = A.coeff(0, 0) + B.coeff(1, 1);
73 coeffMatrix.coeffRef(2, 2) = A.coeff(1, 1) + B.coeff(0, 0);
74 coeffMatrix.coeffRef(3, 3) = A.coeff(1, 1) + B.coeff(1, 1);
75 coeffMatrix.coeffRef(0, 1) = B.coeff(1, 0);
76 coeffMatrix.coeffRef(0, 2) = A.coeff(0, 1);
77 coeffMatrix.coeffRef(1, 0) = B.coeff(0, 1);
78 coeffMatrix.coeffRef(1, 3) = A.coeff(0, 1);
79 coeffMatrix.coeffRef(2, 0) = A.coeff(1, 0);
80 coeffMatrix.coeffRef(2, 3) = B.coeff(1, 0);
81 coeffMatrix.coeffRef(3, 1) = A.coeff(1, 0);
82 coeffMatrix.coeffRef(3, 2) = B.coeff(0, 1);
84 Matrix<Scalar, 4, 1> rhs;
85 rhs.coeffRef(0) = C.coeff(0, 0);
86 rhs.coeffRef(1) = C.coeff(0, 1);
87 rhs.coeffRef(2) = C.coeff(1, 0);
88 rhs.coeffRef(3) = C.coeff(1, 1);
90 Matrix<Scalar, 4, 1> result;
91 result = coeffMatrix.fullPivLu().solve(rhs);
93 X.coeffRef(0, 0) = result.coeff(0);
94 X.coeffRef(0, 1) = result.coeff(1);
95 X.coeffRef(1, 0) = result.coeff(2);
96 X.coeffRef(1, 1) = result.coeff(3);
100template <
typename MatrixType,
typename ResultType>
101void matrix_sqrt_quasi_triangular_2x2_off_diagonal_block(
const MatrixType& T,
Index i,
Index j, ResultType& sqrtT) {
102 typedef typename traits<MatrixType>::Scalar Scalar;
103 Matrix<Scalar, 2, 2> A = sqrtT.template block<2, 2>(i, i);
104 Matrix<Scalar, 2, 2> B = sqrtT.template block<2, 2>(j, j);
105 Matrix<Scalar, 2, 2> C = T.template block<2, 2>(i, j);
106 if (j - i > 2) C -= sqrtT.block(i, i + 2, 2, j - i - 2) * sqrtT.block(i + 2, j, j - i - 2, 2);
107 Matrix<Scalar, 2, 2> X;
108 matrix_sqrt_quasi_triangular_solve_auxiliary_equation(X, A, B, C);
109 sqrtT.template block<2, 2>(i, j) = X;
114template <
typename MatrixType,
typename ResultType>
115void matrix_sqrt_quasi_triangular_diagonal(
const MatrixType& T, ResultType& sqrtT) {
117 const Index size = T.rows();
118 for (
Index i = 0; i < size; i++) {
119 if (i == size - 1 || T.coeff(i + 1, i) == 0) {
120 eigen_assert(T(i, i) >= 0);
121 sqrtT.coeffRef(i, i) =
sqrt(T.coeff(i, i));
123 matrix_sqrt_quasi_triangular_2x2_diagonal_block(T, i, sqrtT);
131template <
typename MatrixType,
typename ResultType>
132void matrix_sqrt_quasi_triangular_off_diagonal(
const MatrixType& T, ResultType& sqrtT) {
133 const Index size = T.rows();
134 for (
Index j = 1; j < size; j++) {
135 if (T.coeff(j, j - 1) != 0)
137 for (
Index i = j - 1; i >= 0; i--) {
138 if (i > 0 && T.coeff(i, i - 1) != 0)
140 bool iBlockIs2x2 = (i < size - 1) && (T.coeff(i + 1, i) != 0);
141 bool jBlockIs2x2 = (j < size - 1) && (T.coeff(j + 1, j) != 0);
142 if (iBlockIs2x2 && jBlockIs2x2)
143 matrix_sqrt_quasi_triangular_2x2_off_diagonal_block(T, i, j, sqrtT);
144 else if (iBlockIs2x2 && !jBlockIs2x2)
145 matrix_sqrt_quasi_triangular_2x1_off_diagonal_block(T, i, j, sqrtT);
146 else if (!iBlockIs2x2 && jBlockIs2x2)
147 matrix_sqrt_quasi_triangular_1x2_off_diagonal_block(T, i, j, sqrtT);
148 else if (!iBlockIs2x2 && !jBlockIs2x2)
149 matrix_sqrt_quasi_triangular_1x1_off_diagonal_block(T, i, j, sqrtT);
171template <
typename MatrixType,
typename ResultType>
173 eigen_assert(
arg.rows() ==
arg.cols());
174 result.resize(
arg.rows(),
arg.cols());
175 internal::matrix_sqrt_quasi_triangular_diagonal(
arg, result);
176 internal::matrix_sqrt_quasi_triangular_off_diagonal(
arg, result);
193template <
typename MatrixType,
typename ResultType>
196 typedef typename MatrixType::Scalar
Scalar;
198 eigen_assert(
arg.rows() ==
arg.cols());
202 result.resize(
arg.rows(),
arg.cols());
203 for (
Index i = 0; i <
arg.rows(); i++) {
204 result.coeffRef(i, i) =
sqrt(
arg.coeff(i, i));
206 for (
Index j = 1; j <
arg.cols(); j++) {
207 for (
Index i = j - 1; i >= 0; i--) {
209 Scalar tmp = (result.row(i).segment(i + 1, j - i - 1) * result.col(j).segment(i + 1, j - i - 1)).value();
211 result.coeffRef(i, j) = (
arg.coeff(i, j) - tmp) / (result.coeff(i, i) + result.coeff(j, j));
225template <typename MatrixType, int IsComplex = NumTraits<typename internal::traits<MatrixType>::Scalar>::IsComplex>
234 template <
typename ResultType>
235 static void run(
const MatrixType&
arg, ResultType& result);
240template <
typename MatrixType>
242 typedef typename MatrixType::PlainObject PlainType;
243 template <
typename ResultType>
244 static void run(
const MatrixType&
arg, ResultType& result) {
245 eigen_assert(
arg.rows() ==
arg.cols());
249 const PlainType& T = schurOfA.matrixT();
250 const PlainType& U = schurOfA.matrixU();
253 PlainType sqrtT = PlainType::Zero(
arg.rows(),
arg.cols());
257 result = U * sqrtT * U.adjoint();
263template <
typename MatrixType>
264struct matrix_sqrt_compute<MatrixType, 1> {
265 typedef typename MatrixType::PlainObject PlainType;
266 template <
typename ResultType>
267 static void run(
const MatrixType&
arg, ResultType& result) {
268 eigen_assert(
arg.rows() ==
arg.cols());
271 const ComplexSchur<PlainType> schurOfA(
arg);
272 const PlainType& T = schurOfA.matrixT();
273 const PlainType& U = schurOfA.matrixU();
280 result = U * (sqrtT.template triangularView<Upper>() * U.adjoint());
298template <
typename Derived>
301 typedef typename internal::ref_selector<Derived>::type DerivedNested;
316 template <
typename ResultType>
317 inline void evalTo(ResultType& result)
const {
318 typedef typename internal::nested_eval<Derived, 10>::type DerivedEvalType;
319 typedef internal::remove_all_t<DerivedEvalType> DerivedEvalTypeClean;
320 DerivedEvalType tmp(m_src);
324 Index rows()
const {
return m_src.rows(); }
325 Index cols()
const {
return m_src.cols(); }
328 const DerivedNested m_src;
332template <
typename Derived>
333struct traits<MatrixSquareRootReturnValue<Derived> > {
334 typedef typename Derived::PlainObject ReturnType;
338template <
typename Derived>
340 eigen_assert(rows() == cols());
const MatrixSquareRootReturnValue< Derived > sqrt() const
Definition MatrixSquareRoot.h:339
Proxy for the matrix square root of some matrix (expression).
Definition MatrixSquareRoot.h:299
void evalTo(ResultType &result) const
Compute the matrix square root.
Definition MatrixSquareRoot.h:317
MatrixSquareRootReturnValue(const Derived &src)
Constructor.
Definition MatrixSquareRoot.h:309
void matrix_sqrt_quasi_triangular(const MatrixType &arg, ResultType &result)
Compute matrix square root of quasi-triangular matrix.
Definition MatrixSquareRoot.h:172
void matrix_sqrt_triangular(const MatrixType &arg, ResultType &result)
Compute matrix square root of triangular matrix.
Definition MatrixSquareRoot.h:194
Namespace containing all symbols from the Eigen library.
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_arg_op< typename Derived::Scalar >, const Derived > arg(const Eigen::ArrayBase< Derived > &x)
Helper struct for computing matrix square roots of general matrices.
Definition MatrixSquareRoot.h:226
static void run(const MatrixType &arg, ResultType &result)
Compute the matrix square root.