MatrixSquareRoot.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2011, 2013 Jitse Niesen <jitse@maths.leeds.ac.uk>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_MATRIX_SQUARE_ROOT
11 #define EIGEN_MATRIX_SQUARE_ROOT
12 
13 // IWYU pragma: private
14 #include "./InternalHeaderCheck.h"
15 
16 namespace Eigen {
17 
18 namespace internal {
19 
20 // pre: T.block(i,i,2,2) has complex conjugate eigenvalues
21 // post: sqrtT.block(i,i,2,2) is square root of T.block(i,i,2,2)
22 template <typename MatrixType, typename ResultType>
24  // TODO: This case (2-by-2 blocks with complex conjugate eigenvalues) is probably hidden somewhere
25  // in EigenSolver. If we expose it, we could call it directly from here.
26  typedef typename traits<MatrixType>::Scalar Scalar;
27  Matrix<Scalar, 2, 2> block = T.template block<2, 2>(i, i);
29  sqrtT.template block<2, 2>(i, i) =
30  (es.eigenvectors() * es.eigenvalues().cwiseSqrt().asDiagonal() * es.eigenvectors().inverse()).real();
31 }
32 
33 // pre: block structure of T is such that (i,j) is a 1x1 block,
34 // all blocks of sqrtT to left of and below (i,j) are correct
35 // post: sqrtT(i,j) has the correct value
36 template <typename MatrixType, typename ResultType>
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));
41 }
42 
43 // similar to compute1x1offDiagonalBlock()
44 template <typename MatrixType, typename ResultType>
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);
50  A += sqrtT.template block<2, 2>(j, j).transpose();
51  sqrtT.template block<1, 2>(i, j).transpose() = A.fullPivLu().solve(rhs.transpose());
52 }
53 
54 // similar to compute1x1offDiagonalBlock()
55 template <typename MatrixType, typename ResultType>
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);
61  A += sqrtT.template block<2, 2>(i, i);
62  sqrtT.template block<2, 1>(i, j) = A.fullPivLu().solve(rhs);
63 }
64 
65 // solves the equation A X + X B = C where all matrices are 2-by-2
66 template <typename MatrixType>
68  const MatrixType& C) {
69  typedef typename traits<MatrixType>::Scalar Scalar;
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);
83 
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);
89 
90  Matrix<Scalar, 4, 1> result;
91  result = coeffMatrix.fullPivLu().solve(rhs);
92 
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);
97 }
98 
99 // similar to compute1x1offDiagonalBlock()
100 template <typename MatrixType, typename ResultType>
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);
109  sqrtT.template block<2, 2>(i, j) = X;
110 }
111 
112 // pre: T is quasi-upper-triangular and sqrtT is a zero matrix of the same size
113 // post: the diagonal blocks of sqrtT are the square roots of the diagonal blocks of T
114 template <typename MatrixType, typename ResultType>
115 void matrix_sqrt_quasi_triangular_diagonal(const MatrixType& T, ResultType& sqrtT) {
116  using std::sqrt;
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));
122  } else {
124  ++i;
125  }
126  }
127 }
128 
129 // pre: T is quasi-upper-triangular and diagonal blocks of sqrtT are square root of diagonal blocks of T.
130 // post: sqrtT is the square root of T.
131 template <typename MatrixType, typename ResultType>
132 void 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) // if T(j-1:j, j-1:j) is a 2-by-2 block
136  continue;
137  for (Index i = j - 1; i >= 0; i--) {
138  if (i > 0 && T.coeff(i, i - 1) != 0) // if T(i-1:i, i-1:i) is a 2-by-2 block
139  continue;
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)
144  else if (iBlockIs2x2 && !jBlockIs2x2)
146  else if (!iBlockIs2x2 && jBlockIs2x2)
148  else if (!iBlockIs2x2 && !jBlockIs2x2)
150  }
151  }
152 }
153 
154 } // end of namespace internal
155 
171 template <typename MatrixType, typename ResultType>
172 void matrix_sqrt_quasi_triangular(const MatrixType& arg, ResultType& result) {
173  eigen_assert(arg.rows() == arg.cols());
174  result.resize(arg.rows(), arg.cols());
177 }
178 
193 template <typename MatrixType, typename ResultType>
194 void matrix_sqrt_triangular(const MatrixType& arg, ResultType& result) {
195  using std::sqrt;
196  typedef typename MatrixType::Scalar Scalar;
197 
198  eigen_assert(arg.rows() == arg.cols());
199 
200  // Compute square root of arg and store it in upper triangular part of result
201  // This uses that the square root of triangular matrices can be computed directly.
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));
205  }
206  for (Index j = 1; j < arg.cols(); j++) {
207  for (Index i = j - 1; i >= 0; i--) {
208  // if i = j-1, then segment has length 0 so tmp = 0
209  Scalar tmp = (result.row(i).segment(i + 1, j - i - 1) * result.col(j).segment(i + 1, j - i - 1)).value();
210  // denominator may be zero if original matrix is singular
211  result.coeffRef(i, j) = (arg.coeff(i, j) - tmp) / (result.coeff(i, i) + result.coeff(j, j));
212  }
213  }
214 }
215 
216 namespace internal {
217 
234  template <typename ResultType>
235  static void run(const MatrixType& arg, ResultType& result);
236 };
237 
238 // ********** Partial specialization for real matrices **********
239 
240 template <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());
246 
247  // Compute Schur decomposition of arg
248  const RealSchur<PlainType> schurOfA(arg);
249  const PlainType& T = schurOfA.matrixT();
250  const PlainType& U = schurOfA.matrixU();
251 
252  // Compute square root of T
253  PlainType sqrtT = PlainType::Zero(arg.rows(), arg.cols());
255 
256  // Compute square root of arg
257  result = U * sqrtT * U.adjoint();
258  }
259 };
260 
261 // ********** Partial specialization for complex matrices **********
262 
263 template <typename MatrixType>
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());
269 
270  // Compute Schur decomposition of arg
272  const PlainType& T = schurOfA.matrixT();
273  const PlainType& U = schurOfA.matrixU();
274 
275  // Compute square root of T
276  PlainType sqrtT;
277  matrix_sqrt_triangular(T, sqrtT);
278 
279  // Compute square root of arg
280  result = U * (sqrtT.template triangularView<Upper>() * U.adjoint());
281  }
282 };
283 
284 } // end namespace internal
285 
298 template <typename Derived>
299 class MatrixSquareRootReturnValue : public ReturnByValue<MatrixSquareRootReturnValue<Derived> > {
300  protected:
302 
303  public:
309  explicit MatrixSquareRootReturnValue(const Derived& src) : m_src(src) {}
310 
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);
322  }
323 
324  Index rows() const { return m_src.rows(); }
325  Index cols() const { return m_src.cols(); }
326 
327  protected:
329 };
330 
331 namespace internal {
332 template <typename Derived>
334  typedef typename Derived::PlainObject ReturnType;
335 };
336 } // namespace internal
337 
338 template <typename Derived>
340  eigen_assert(rows() == cols());
341  return MatrixSquareRootReturnValue<Derived>(derived());
342 }
343 
344 } // end namespace Eigen
345 
346 #endif // EIGEN_MATRIX_FUNCTION
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
int i
Definition: BiCGSTAB_step_by_step.cpp:9
cout<< "Here is a random 4x4 matrix, A:"<< endl<< A<< endl<< endl;ComplexSchur< MatrixXcf > schurOfA(A, false)
EigenSolver< MatrixXf > es
Definition: EigenSolver_compute.cpp:1
Eigen::Triplet< double > T
Definition: EigenUnitTest.cpp:11
#define eigen_assert(x)
Definition: Macros.h:910
m m block(1, 0, 2, 2)<< 4
int rows
Definition: Tutorial_commainit_02.cpp:1
int cols
Definition: Tutorial_commainit_02.cpp:1
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
Performs a complex Schur decomposition of a real or complex square matrix.
Definition: ComplexSchur.h:56
Computes eigenvalues and eigenvectors of general matrices.
Definition: EigenSolver.h:68
Proxy for the matrix square root of some matrix (expression).
Definition: MatrixSquareRoot.h:299
internal::ref_selector< Derived >::type DerivedNested
Definition: MatrixSquareRoot.h:301
Index cols() const
Definition: MatrixSquareRoot.h:325
const DerivedNested m_src
Definition: MatrixSquareRoot.h:328
void evalTo(ResultType &result) const
Compute the matrix square root.
Definition: MatrixSquareRoot.h:317
MatrixSquareRootReturnValue(const Derived &src)
Constructor.
Definition: MatrixSquareRoot.h:309
Index rows() const
Definition: MatrixSquareRoot.h:324
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE Scalar & coeffRef(Index rowId, Index colId)
Definition: PlainObjectBase.h:217
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE const Scalar & coeff(Index rowId, Index colId) const
Definition: PlainObjectBase.h:198
Performs a real Schur decomposition of a square matrix.
Definition: RealSchur.h:58
Definition: ReturnByValue.h:50
Definition: matrices.h:74
@ IsComplex
Definition: common.h:73
#define X
Definition: icosphere.cpp:20
Eigen::Matrix< Scalar, Dynamic, Dynamic, ColMajor > tmp
Definition: level3_impl.h:365
void matrix_sqrt_quasi_triangular_1x2_off_diagonal_block(const MatrixType &T, Index i, Index j, ResultType &sqrtT)
Definition: MatrixSquareRoot.h:45
void matrix_sqrt_quasi_triangular_solve_auxiliary_equation(MatrixType &X, const MatrixType &A, const MatrixType &B, const MatrixType &C)
Definition: MatrixSquareRoot.h:67
void matrix_sqrt_quasi_triangular_diagonal(const MatrixType &T, ResultType &sqrtT)
Definition: MatrixSquareRoot.h:115
void matrix_sqrt_quasi_triangular_2x2_diagonal_block(const MatrixType &T, Index i, ResultType &sqrtT)
Definition: MatrixSquareRoot.h:23
typename remove_all< T >::type remove_all_t
Definition: Meta.h:142
void matrix_sqrt_quasi_triangular_off_diagonal(const MatrixType &T, ResultType &sqrtT)
Definition: MatrixSquareRoot.h:132
void matrix_sqrt_quasi_triangular_2x2_off_diagonal_block(const MatrixType &T, Index i, Index j, ResultType &sqrtT)
Definition: MatrixSquareRoot.h:101
void matrix_sqrt_quasi_triangular_1x1_off_diagonal_block(const MatrixType &T, Index i, Index j, ResultType &sqrtT)
Definition: MatrixSquareRoot.h:37
void matrix_sqrt_quasi_triangular_2x1_off_diagonal_block(const MatrixType &T, Index i, Index j, ResultType &sqrtT)
Definition: MatrixSquareRoot.h:56
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
squared absolute value
Definition: GlobalFunctions.h:87
void matrix_sqrt_quasi_triangular(const MatrixType &arg, ResultType &result)
Compute matrix square root of quasi-triangular matrix.
Definition: MatrixSquareRoot.h:172
const AutoDiffScalar< DerType > & real(const AutoDiffScalar< DerType > &x)
Definition: AutoDiffScalar.h:486
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
void matrix_sqrt_triangular(const MatrixType &arg, ResultType &result)
Compute matrix square root of triangular matrix.
Definition: MatrixSquareRoot.h:194
double U
Swimming speed.
Definition: two_d_variable_diff_adapt.cc:53
Definition: Eigen_Colamd.h:49
double Zero
Definition: pseudosolid_node_update_elements.cc:35
MatrixType::PlainObject PlainType
Definition: MatrixSquareRoot.h:242
static void run(const MatrixType &arg, ResultType &result)
Definition: MatrixSquareRoot.h:244
static void run(const MatrixType &arg, ResultType &result)
Definition: MatrixSquareRoot.h:267
MatrixType::PlainObject PlainType
Definition: MatrixSquareRoot.h:265
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.
std::conditional_t< Evaluate, PlainObject, typename ref_selector< T >::type > type
Definition: XprHelper.h:549
std::conditional_t< bool(traits< T >::Flags &NestByRefBit), T const &, const T > type
Definition: XprHelper.h:507
Derived::PlainObject ReturnType
Definition: MatrixSquareRoot.h:334
Definition: ForwardDeclarations.h:21
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2