Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 1 | // This file is part of Eigen, a lightweight C++ template library |
| 2 | // for linear algebra. |
| 3 | // |
| 4 | // Copyright (C) 2009 Jitse Niesen <jitse@maths.leeds.ac.uk> |
| 5 | // Copyright (C) 2012 Chen-Pang He <jdh8@ms63.hinet.net> |
| 6 | // |
| 7 | // This Source Code Form is subject to the terms of the Mozilla |
| 8 | // Public License v. 2.0. If a copy of the MPL was not distributed |
| 9 | // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. |
| 10 | |
| 11 | #ifndef EIGEN_MATRIX_FUNCTIONS |
| 12 | #define EIGEN_MATRIX_FUNCTIONS |
| 13 | |
| 14 | #include <cfloat> |
| 15 | #include <list> |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 16 | |
| 17 | #include <Eigen/Core> |
| 18 | #include <Eigen/LU> |
| 19 | #include <Eigen/Eigenvalues> |
| 20 | |
| 21 | /** |
| 22 | * \defgroup MatrixFunctions_Module Matrix functions module |
| 23 | * \brief This module aims to provide various methods for the computation of |
| 24 | * matrix functions. |
| 25 | * |
| 26 | * To use this module, add |
| 27 | * \code |
| 28 | * #include <unsupported/Eigen/MatrixFunctions> |
| 29 | * \endcode |
| 30 | * at the start of your source file. |
| 31 | * |
| 32 | * This module defines the following MatrixBase methods. |
| 33 | * - \ref matrixbase_cos "MatrixBase::cos()", for computing the matrix cosine |
| 34 | * - \ref matrixbase_cosh "MatrixBase::cosh()", for computing the matrix hyperbolic cosine |
| 35 | * - \ref matrixbase_exp "MatrixBase::exp()", for computing the matrix exponential |
| 36 | * - \ref matrixbase_log "MatrixBase::log()", for computing the matrix logarithm |
| 37 | * - \ref matrixbase_pow "MatrixBase::pow()", for computing the matrix power |
| 38 | * - \ref matrixbase_matrixfunction "MatrixBase::matrixFunction()", for computing general matrix functions |
| 39 | * - \ref matrixbase_sin "MatrixBase::sin()", for computing the matrix sine |
| 40 | * - \ref matrixbase_sinh "MatrixBase::sinh()", for computing the matrix hyperbolic sine |
| 41 | * - \ref matrixbase_sqrt "MatrixBase::sqrt()", for computing the matrix square root |
| 42 | * |
| 43 | * These methods are the main entry points to this module. |
| 44 | * |
| 45 | * %Matrix functions are defined as follows. Suppose that \f$ f \f$ |
| 46 | * is an entire function (that is, a function on the complex plane |
| 47 | * that is everywhere complex differentiable). Then its Taylor |
| 48 | * series |
| 49 | * \f[ f(0) + f'(0) x + \frac{f''(0)}{2} x^2 + \frac{f'''(0)}{3!} x^3 + \cdots \f] |
| 50 | * converges to \f$ f(x) \f$. In this case, we can define the matrix |
| 51 | * function by the same series: |
| 52 | * \f[ f(M) = f(0) + f'(0) M + \frac{f''(0)}{2} M^2 + \frac{f'''(0)}{3!} M^3 + \cdots \f] |
| 53 | * |
| 54 | */ |
| 55 | |
| 56 | #include "src/MatrixFunctions/MatrixExponential.h" |
| 57 | #include "src/MatrixFunctions/MatrixFunction.h" |
| 58 | #include "src/MatrixFunctions/MatrixSquareRoot.h" |
| 59 | #include "src/MatrixFunctions/MatrixLogarithm.h" |
| 60 | #include "src/MatrixFunctions/MatrixPower.h" |
| 61 | |
| 62 | |
| 63 | /** |
| 64 | \page matrixbaseextra_page |
| 65 | \ingroup MatrixFunctions_Module |
| 66 | |
| 67 | \section matrixbaseextra MatrixBase methods defined in the MatrixFunctions module |
| 68 | |
| 69 | The remainder of the page documents the following MatrixBase methods |
| 70 | which are defined in the MatrixFunctions module. |
| 71 | |
| 72 | |
| 73 | |
| 74 | \subsection matrixbase_cos MatrixBase::cos() |
| 75 | |
| 76 | Compute the matrix cosine. |
| 77 | |
| 78 | \code |
| 79 | const MatrixFunctionReturnValue<Derived> MatrixBase<Derived>::cos() const |
| 80 | \endcode |
| 81 | |
| 82 | \param[in] M a square matrix. |
| 83 | \returns expression representing \f$ \cos(M) \f$. |
| 84 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 85 | This function computes the matrix cosine. Use ArrayBase::cos() for computing the entry-wise cosine. |
| 86 | |
| 87 | The implementation calls \ref matrixbase_matrixfunction "matrixFunction()" with StdStemFunctions::cos(). |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 88 | |
| 89 | \sa \ref matrixbase_sin "sin()" for an example. |
| 90 | |
| 91 | |
| 92 | |
| 93 | \subsection matrixbase_cosh MatrixBase::cosh() |
| 94 | |
| 95 | Compute the matrix hyberbolic cosine. |
| 96 | |
| 97 | \code |
| 98 | const MatrixFunctionReturnValue<Derived> MatrixBase<Derived>::cosh() const |
| 99 | \endcode |
| 100 | |
| 101 | \param[in] M a square matrix. |
| 102 | \returns expression representing \f$ \cosh(M) \f$ |
| 103 | |
| 104 | This function calls \ref matrixbase_matrixfunction "matrixFunction()" with StdStemFunctions::cosh(). |
| 105 | |
| 106 | \sa \ref matrixbase_sinh "sinh()" for an example. |
| 107 | |
| 108 | |
| 109 | |
| 110 | \subsection matrixbase_exp MatrixBase::exp() |
| 111 | |
| 112 | Compute the matrix exponential. |
| 113 | |
| 114 | \code |
| 115 | const MatrixExponentialReturnValue<Derived> MatrixBase<Derived>::exp() const |
| 116 | \endcode |
| 117 | |
| 118 | \param[in] M matrix whose exponential is to be computed. |
| 119 | \returns expression representing the matrix exponential of \p M. |
| 120 | |
| 121 | The matrix exponential of \f$ M \f$ is defined by |
| 122 | \f[ \exp(M) = \sum_{k=0}^\infty \frac{M^k}{k!}. \f] |
| 123 | The matrix exponential can be used to solve linear ordinary |
| 124 | differential equations: the solution of \f$ y' = My \f$ with the |
| 125 | initial condition \f$ y(0) = y_0 \f$ is given by |
| 126 | \f$ y(t) = \exp(M) y_0 \f$. |
| 127 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 128 | The matrix exponential is different from applying the exp function to all the entries in the matrix. |
| 129 | Use ArrayBase::exp() if you want to do the latter. |
| 130 | |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 131 | The cost of the computation is approximately \f$ 20 n^3 \f$ for |
| 132 | matrices of size \f$ n \f$. The number 20 depends weakly on the |
| 133 | norm of the matrix. |
| 134 | |
| 135 | The matrix exponential is computed using the scaling-and-squaring |
| 136 | method combined with Padé approximation. The matrix is first |
| 137 | rescaled, then the exponential of the reduced matrix is computed |
| 138 | approximant, and then the rescaling is undone by repeated |
| 139 | squaring. The degree of the Padé approximant is chosen such |
| 140 | that the approximation error is less than the round-off |
| 141 | error. However, errors may accumulate during the squaring phase. |
| 142 | |
| 143 | Details of the algorithm can be found in: Nicholas J. Higham, "The |
| 144 | scaling and squaring method for the matrix exponential revisited," |
| 145 | <em>SIAM J. %Matrix Anal. Applic.</em>, <b>26</b>:1179–1193, |
| 146 | 2005. |
| 147 | |
| 148 | Example: The following program checks that |
| 149 | \f[ \exp \left[ \begin{array}{ccc} |
| 150 | 0 & \frac14\pi & 0 \\ |
| 151 | -\frac14\pi & 0 & 0 \\ |
| 152 | 0 & 0 & 0 |
| 153 | \end{array} \right] = \left[ \begin{array}{ccc} |
| 154 | \frac12\sqrt2 & -\frac12\sqrt2 & 0 \\ |
| 155 | \frac12\sqrt2 & \frac12\sqrt2 & 0 \\ |
| 156 | 0 & 0 & 1 |
| 157 | \end{array} \right]. \f] |
| 158 | This corresponds to a rotation of \f$ \frac14\pi \f$ radians around |
| 159 | the z-axis. |
| 160 | |
| 161 | \include MatrixExponential.cpp |
| 162 | Output: \verbinclude MatrixExponential.out |
| 163 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 164 | \note \p M has to be a matrix of \c float, \c double, `long double` |
| 165 | \c complex<float>, \c complex<double>, or `complex<long double>` . |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 166 | |
| 167 | |
| 168 | \subsection matrixbase_log MatrixBase::log() |
| 169 | |
| 170 | Compute the matrix logarithm. |
| 171 | |
| 172 | \code |
| 173 | const MatrixLogarithmReturnValue<Derived> MatrixBase<Derived>::log() const |
| 174 | \endcode |
| 175 | |
| 176 | \param[in] M invertible matrix whose logarithm is to be computed. |
| 177 | \returns expression representing the matrix logarithm root of \p M. |
| 178 | |
| 179 | The matrix logarithm of \f$ M \f$ is a matrix \f$ X \f$ such that |
| 180 | \f$ \exp(X) = M \f$ where exp denotes the matrix exponential. As for |
| 181 | the scalar logarithm, the equation \f$ \exp(X) = M \f$ may have |
| 182 | multiple solutions; this function returns a matrix whose eigenvalues |
| 183 | have imaginary part in the interval \f$ (-\pi,\pi] \f$. |
| 184 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 185 | The matrix logarithm is different from applying the log function to all the entries in the matrix. |
| 186 | Use ArrayBase::log() if you want to do the latter. |
| 187 | |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 188 | In the real case, the matrix \f$ M \f$ should be invertible and |
| 189 | it should have no eigenvalues which are real and negative (pairs of |
| 190 | complex conjugate eigenvalues are allowed). In the complex case, it |
| 191 | only needs to be invertible. |
| 192 | |
| 193 | This function computes the matrix logarithm using the Schur-Parlett |
| 194 | algorithm as implemented by MatrixBase::matrixFunction(). The |
| 195 | logarithm of an atomic block is computed by MatrixLogarithmAtomic, |
| 196 | which uses direct computation for 1-by-1 and 2-by-2 blocks and an |
| 197 | inverse scaling-and-squaring algorithm for bigger blocks, with the |
| 198 | square roots computed by MatrixBase::sqrt(). |
| 199 | |
| 200 | Details of the algorithm can be found in Section 11.6.2 of: |
| 201 | Nicholas J. Higham, |
| 202 | <em>Functions of Matrices: Theory and Computation</em>, |
| 203 | SIAM 2008. ISBN 978-0-898716-46-7. |
| 204 | |
| 205 | Example: The following program checks that |
| 206 | \f[ \log \left[ \begin{array}{ccc} |
| 207 | \frac12\sqrt2 & -\frac12\sqrt2 & 0 \\ |
| 208 | \frac12\sqrt2 & \frac12\sqrt2 & 0 \\ |
| 209 | 0 & 0 & 1 |
| 210 | \end{array} \right] = \left[ \begin{array}{ccc} |
| 211 | 0 & \frac14\pi & 0 \\ |
| 212 | -\frac14\pi & 0 & 0 \\ |
| 213 | 0 & 0 & 0 |
| 214 | \end{array} \right]. \f] |
| 215 | This corresponds to a rotation of \f$ \frac14\pi \f$ radians around |
| 216 | the z-axis. This is the inverse of the example used in the |
| 217 | documentation of \ref matrixbase_exp "exp()". |
| 218 | |
| 219 | \include MatrixLogarithm.cpp |
| 220 | Output: \verbinclude MatrixLogarithm.out |
| 221 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 222 | \note \p M has to be a matrix of \c float, \c double, `long |
| 223 | double`, \c complex<float>, \c complex<double>, or `complex<long double>`. |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 224 | |
| 225 | \sa MatrixBase::exp(), MatrixBase::matrixFunction(), |
| 226 | class MatrixLogarithmAtomic, MatrixBase::sqrt(). |
| 227 | |
| 228 | |
| 229 | \subsection matrixbase_pow MatrixBase::pow() |
| 230 | |
| 231 | Compute the matrix raised to arbitrary real power. |
| 232 | |
| 233 | \code |
| 234 | const MatrixPowerReturnValue<Derived> MatrixBase<Derived>::pow(RealScalar p) const |
| 235 | \endcode |
| 236 | |
| 237 | \param[in] M base of the matrix power, should be a square matrix. |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 238 | \param[in] p exponent of the matrix power. |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 239 | |
| 240 | The matrix power \f$ M^p \f$ is defined as \f$ \exp(p \log(M)) \f$, |
| 241 | where exp denotes the matrix exponential, and log denotes the matrix |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 242 | logarithm. This is different from raising all the entries in the matrix |
| 243 | to the p-th power. Use ArrayBase::pow() if you want to do the latter. |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 244 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 245 | If \p p is complex, the scalar type of \p M should be the type of \p |
| 246 | p . \f$ M^p \f$ simply evaluates into \f$ \exp(p \log(M)) \f$. |
| 247 | Therefore, the matrix \f$ M \f$ should meet the conditions to be an |
| 248 | argument of matrix logarithm. |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 249 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 250 | If \p p is real, it is casted into the real scalar type of \p M. Then |
| 251 | this function computes the matrix power using the Schur-Padé |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 252 | algorithm as implemented by class MatrixPower. The exponent is split |
| 253 | into integral part and fractional part, where the fractional part is |
| 254 | in the interval \f$ (-1, 1) \f$. The main diagonal and the first |
| 255 | super-diagonal is directly computed. |
| 256 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 257 | If \p M is singular with a semisimple zero eigenvalue and \p p is |
| 258 | positive, the Schur factor \f$ T \f$ is reordered with Givens |
| 259 | rotations, i.e. |
| 260 | |
| 261 | \f[ T = \left[ \begin{array}{cc} |
| 262 | T_1 & T_2 \\ |
| 263 | 0 & 0 |
| 264 | \end{array} \right] \f] |
| 265 | |
| 266 | where \f$ T_1 \f$ is invertible. Then \f$ T^p \f$ is given by |
| 267 | |
| 268 | \f[ T^p = \left[ \begin{array}{cc} |
| 269 | T_1^p & T_1^{-1} T_1^p T_2 \\ |
| 270 | 0 & 0 |
| 271 | \end{array}. \right] \f] |
| 272 | |
| 273 | \warning Fractional power of a matrix with a non-semisimple zero |
| 274 | eigenvalue is not well-defined. We introduce an assertion failure |
| 275 | against inaccurate result, e.g. \code |
| 276 | #include <unsupported/Eigen/MatrixFunctions> |
| 277 | #include <iostream> |
| 278 | |
| 279 | int main() |
| 280 | { |
| 281 | Eigen::Matrix4d A; |
| 282 | A << 0, 0, 2, 3, |
| 283 | 0, 0, 4, 5, |
| 284 | 0, 0, 6, 7, |
| 285 | 0, 0, 8, 9; |
| 286 | std::cout << A.pow(0.37) << std::endl; |
| 287 | |
| 288 | // The 1 makes eigenvalue 0 non-semisimple. |
| 289 | A.coeffRef(0, 1) = 1; |
| 290 | |
| 291 | // This fails if EIGEN_NO_DEBUG is undefined. |
| 292 | std::cout << A.pow(0.37) << std::endl; |
| 293 | |
| 294 | return 0; |
| 295 | } |
| 296 | \endcode |
| 297 | |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 298 | Details of the algorithm can be found in: Nicholas J. Higham and |
| 299 | Lijing Lin, "A Schur-Padé algorithm for fractional powers of a |
| 300 | matrix," <em>SIAM J. %Matrix Anal. Applic.</em>, |
| 301 | <b>32(3)</b>:1056–1078, 2011. |
| 302 | |
| 303 | Example: The following program checks that |
| 304 | \f[ \left[ \begin{array}{ccc} |
| 305 | \cos1 & -\sin1 & 0 \\ |
| 306 | \sin1 & \cos1 & 0 \\ |
| 307 | 0 & 0 & 1 |
| 308 | \end{array} \right]^{\frac14\pi} = \left[ \begin{array}{ccc} |
| 309 | \frac12\sqrt2 & -\frac12\sqrt2 & 0 \\ |
| 310 | \frac12\sqrt2 & \frac12\sqrt2 & 0 \\ |
| 311 | 0 & 0 & 1 |
| 312 | \end{array} \right]. \f] |
| 313 | This corresponds to \f$ \frac14\pi \f$ rotations of 1 radian around |
| 314 | the z-axis. |
| 315 | |
| 316 | \include MatrixPower.cpp |
| 317 | Output: \verbinclude MatrixPower.out |
| 318 | |
| 319 | MatrixBase::pow() is user-friendly. However, there are some |
| 320 | circumstances under which you should use class MatrixPower directly. |
| 321 | MatrixPower can save the result of Schur decomposition, so it's |
| 322 | better for computing various powers for the same matrix. |
| 323 | |
| 324 | Example: |
| 325 | \include MatrixPower_optimal.cpp |
| 326 | Output: \verbinclude MatrixPower_optimal.out |
| 327 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 328 | \note \p M has to be a matrix of \c float, \c double, `long |
| 329 | double`, \c complex<float>, \c complex<double>, or |
| 330 | \c complex<long double> . |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 331 | |
| 332 | \sa MatrixBase::exp(), MatrixBase::log(), class MatrixPower. |
| 333 | |
| 334 | |
| 335 | \subsection matrixbase_matrixfunction MatrixBase::matrixFunction() |
| 336 | |
| 337 | Compute a matrix function. |
| 338 | |
| 339 | \code |
| 340 | const MatrixFunctionReturnValue<Derived> MatrixBase<Derived>::matrixFunction(typename internal::stem_function<typename internal::traits<Derived>::Scalar>::type f) const |
| 341 | \endcode |
| 342 | |
| 343 | \param[in] M argument of matrix function, should be a square matrix. |
| 344 | \param[in] f an entire function; \c f(x,n) should compute the n-th |
| 345 | derivative of f at x. |
| 346 | \returns expression representing \p f applied to \p M. |
| 347 | |
| 348 | Suppose that \p M is a matrix whose entries have type \c Scalar. |
| 349 | Then, the second argument, \p f, should be a function with prototype |
| 350 | \code |
| 351 | ComplexScalar f(ComplexScalar, int) |
| 352 | \endcode |
| 353 | where \c ComplexScalar = \c std::complex<Scalar> if \c Scalar is |
| 354 | real (e.g., \c float or \c double) and \c ComplexScalar = |
| 355 | \c Scalar if \c Scalar is complex. The return value of \c f(x,n) |
| 356 | should be \f$ f^{(n)}(x) \f$, the n-th derivative of f at x. |
| 357 | |
| 358 | This routine uses the algorithm described in: |
| 359 | Philip Davies and Nicholas J. Higham, |
| 360 | "A Schur-Parlett algorithm for computing matrix functions", |
| 361 | <em>SIAM J. %Matrix Anal. Applic.</em>, <b>25</b>:464–485, 2003. |
| 362 | |
| 363 | The actual work is done by the MatrixFunction class. |
| 364 | |
| 365 | Example: The following program checks that |
| 366 | \f[ \exp \left[ \begin{array}{ccc} |
| 367 | 0 & \frac14\pi & 0 \\ |
| 368 | -\frac14\pi & 0 & 0 \\ |
| 369 | 0 & 0 & 0 |
| 370 | \end{array} \right] = \left[ \begin{array}{ccc} |
| 371 | \frac12\sqrt2 & -\frac12\sqrt2 & 0 \\ |
| 372 | \frac12\sqrt2 & \frac12\sqrt2 & 0 \\ |
| 373 | 0 & 0 & 1 |
| 374 | \end{array} \right]. \f] |
| 375 | This corresponds to a rotation of \f$ \frac14\pi \f$ radians around |
| 376 | the z-axis. This is the same example as used in the documentation |
| 377 | of \ref matrixbase_exp "exp()". |
| 378 | |
| 379 | \include MatrixFunction.cpp |
| 380 | Output: \verbinclude MatrixFunction.out |
| 381 | |
| 382 | Note that the function \c expfn is defined for complex numbers |
| 383 | \c x, even though the matrix \c A is over the reals. Instead of |
| 384 | \c expfn, we could also have used StdStemFunctions::exp: |
| 385 | \code |
| 386 | A.matrixFunction(StdStemFunctions<std::complex<double> >::exp, &B); |
| 387 | \endcode |
| 388 | |
| 389 | |
| 390 | |
| 391 | \subsection matrixbase_sin MatrixBase::sin() |
| 392 | |
| 393 | Compute the matrix sine. |
| 394 | |
| 395 | \code |
| 396 | const MatrixFunctionReturnValue<Derived> MatrixBase<Derived>::sin() const |
| 397 | \endcode |
| 398 | |
| 399 | \param[in] M a square matrix. |
| 400 | \returns expression representing \f$ \sin(M) \f$. |
| 401 | |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 402 | This function computes the matrix sine. Use ArrayBase::sin() for computing the entry-wise sine. |
| 403 | |
| 404 | The implementation calls \ref matrixbase_matrixfunction "matrixFunction()" with StdStemFunctions::sin(). |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 405 | |
| 406 | Example: \include MatrixSine.cpp |
| 407 | Output: \verbinclude MatrixSine.out |
| 408 | |
| 409 | |
| 410 | |
| 411 | \subsection matrixbase_sinh MatrixBase::sinh() |
| 412 | |
| 413 | Compute the matrix hyperbolic sine. |
| 414 | |
| 415 | \code |
| 416 | MatrixFunctionReturnValue<Derived> MatrixBase<Derived>::sinh() const |
| 417 | \endcode |
| 418 | |
| 419 | \param[in] M a square matrix. |
| 420 | \returns expression representing \f$ \sinh(M) \f$ |
| 421 | |
| 422 | This function calls \ref matrixbase_matrixfunction "matrixFunction()" with StdStemFunctions::sinh(). |
| 423 | |
| 424 | Example: \include MatrixSinh.cpp |
| 425 | Output: \verbinclude MatrixSinh.out |
| 426 | |
| 427 | |
| 428 | \subsection matrixbase_sqrt MatrixBase::sqrt() |
| 429 | |
| 430 | Compute the matrix square root. |
| 431 | |
| 432 | \code |
| 433 | const MatrixSquareRootReturnValue<Derived> MatrixBase<Derived>::sqrt() const |
| 434 | \endcode |
| 435 | |
| 436 | \param[in] M invertible matrix whose square root is to be computed. |
| 437 | \returns expression representing the matrix square root of \p M. |
| 438 | |
| 439 | The matrix square root of \f$ M \f$ is the matrix \f$ M^{1/2} \f$ |
| 440 | whose square is the original matrix; so if \f$ S = M^{1/2} \f$ then |
Austin Schuh | 189376f | 2018-12-20 22:11:15 +1100 | [diff] [blame^] | 441 | \f$ S^2 = M \f$. This is different from taking the square root of all |
| 442 | the entries in the matrix; use ArrayBase::sqrt() if you want to do the |
| 443 | latter. |
Brian Silverman | 72890c2 | 2015-09-19 14:37:37 -0400 | [diff] [blame] | 444 | |
| 445 | In the <b>real case</b>, the matrix \f$ M \f$ should be invertible and |
| 446 | it should have no eigenvalues which are real and negative (pairs of |
| 447 | complex conjugate eigenvalues are allowed). In that case, the matrix |
| 448 | has a square root which is also real, and this is the square root |
| 449 | computed by this function. |
| 450 | |
| 451 | The matrix square root is computed by first reducing the matrix to |
| 452 | quasi-triangular form with the real Schur decomposition. The square |
| 453 | root of the quasi-triangular matrix can then be computed directly. The |
| 454 | cost is approximately \f$ 25 n^3 \f$ real flops for the real Schur |
| 455 | decomposition and \f$ 3\frac13 n^3 \f$ real flops for the remainder |
| 456 | (though the computation time in practice is likely more than this |
| 457 | indicates). |
| 458 | |
| 459 | Details of the algorithm can be found in: Nicholas J. Highan, |
| 460 | "Computing real square roots of a real matrix", <em>Linear Algebra |
| 461 | Appl.</em>, 88/89:405–430, 1987. |
| 462 | |
| 463 | If the matrix is <b>positive-definite symmetric</b>, then the square |
| 464 | root is also positive-definite symmetric. In this case, it is best to |
| 465 | use SelfAdjointEigenSolver::operatorSqrt() to compute it. |
| 466 | |
| 467 | In the <b>complex case</b>, the matrix \f$ M \f$ should be invertible; |
| 468 | this is a restriction of the algorithm. The square root computed by |
| 469 | this algorithm is the one whose eigenvalues have an argument in the |
| 470 | interval \f$ (-\frac12\pi, \frac12\pi] \f$. This is the usual branch |
| 471 | cut. |
| 472 | |
| 473 | The computation is the same as in the real case, except that the |
| 474 | complex Schur decomposition is used to reduce the matrix to a |
| 475 | triangular matrix. The theoretical cost is the same. Details are in: |
| 476 | Åke Björck and Sven Hammarling, "A Schur method for the |
| 477 | square root of a matrix", <em>Linear Algebra Appl.</em>, |
| 478 | 52/53:127–140, 1983. |
| 479 | |
| 480 | Example: The following program checks that the square root of |
| 481 | \f[ \left[ \begin{array}{cc} |
| 482 | \cos(\frac13\pi) & -\sin(\frac13\pi) \\ |
| 483 | \sin(\frac13\pi) & \cos(\frac13\pi) |
| 484 | \end{array} \right], \f] |
| 485 | corresponding to a rotation over 60 degrees, is a rotation over 30 degrees: |
| 486 | \f[ \left[ \begin{array}{cc} |
| 487 | \cos(\frac16\pi) & -\sin(\frac16\pi) \\ |
| 488 | \sin(\frac16\pi) & \cos(\frac16\pi) |
| 489 | \end{array} \right]. \f] |
| 490 | |
| 491 | \include MatrixSquareRoot.cpp |
| 492 | Output: \verbinclude MatrixSquareRoot.out |
| 493 | |
| 494 | \sa class RealSchur, class ComplexSchur, class MatrixSquareRoot, |
| 495 | SelfAdjointEigenSolver::operatorSqrt(). |
| 496 | |
| 497 | */ |
| 498 | |
| 499 | #endif // EIGEN_MATRIX_FUNCTIONS |
| 500 | |