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