File indexing completed on 2025-01-18 09:57:05
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010 #ifndef EIGEN_MATRIX_POWER
0011 #define EIGEN_MATRIX_POWER
0012
0013 namespace Eigen {
0014
0015 template<typename MatrixType> class MatrixPower;
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038 template<typename MatrixType>
0039 class MatrixPowerParenthesesReturnValue : public ReturnByValue< MatrixPowerParenthesesReturnValue<MatrixType> >
0040 {
0041 public:
0042 typedef typename MatrixType::RealScalar RealScalar;
0043
0044
0045
0046
0047
0048
0049
0050 MatrixPowerParenthesesReturnValue(MatrixPower<MatrixType>& pow, RealScalar p) : m_pow(pow), m_p(p)
0051 { }
0052
0053
0054
0055
0056
0057
0058 template<typename ResultType>
0059 inline void evalTo(ResultType& result) const
0060 { m_pow.compute(result, m_p); }
0061
0062 Index rows() const { return m_pow.rows(); }
0063 Index cols() const { return m_pow.cols(); }
0064
0065 private:
0066 MatrixPower<MatrixType>& m_pow;
0067 const RealScalar m_p;
0068 };
0069
0070
0071
0072
0073
0074
0075
0076
0077
0078
0079
0080
0081
0082
0083
0084
0085 template<typename MatrixType>
0086 class MatrixPowerAtomic : internal::noncopyable
0087 {
0088 private:
0089 enum {
0090 RowsAtCompileTime = MatrixType::RowsAtCompileTime,
0091 MaxRowsAtCompileTime = MatrixType::MaxRowsAtCompileTime
0092 };
0093 typedef typename MatrixType::Scalar Scalar;
0094 typedef typename MatrixType::RealScalar RealScalar;
0095 typedef std::complex<RealScalar> ComplexScalar;
0096 typedef Block<MatrixType,Dynamic,Dynamic> ResultType;
0097
0098 const MatrixType& m_A;
0099 RealScalar m_p;
0100
0101 void computePade(int degree, const MatrixType& IminusT, ResultType& res) const;
0102 void compute2x2(ResultType& res, RealScalar p) const;
0103 void computeBig(ResultType& res) const;
0104 static int getPadeDegree(float normIminusT);
0105 static int getPadeDegree(double normIminusT);
0106 static int getPadeDegree(long double normIminusT);
0107 static ComplexScalar computeSuperDiag(const ComplexScalar&, const ComplexScalar&, RealScalar p);
0108 static RealScalar computeSuperDiag(RealScalar, RealScalar, RealScalar p);
0109
0110 public:
0111
0112
0113
0114
0115
0116
0117
0118
0119
0120
0121
0122 MatrixPowerAtomic(const MatrixType& T, RealScalar p);
0123
0124
0125
0126
0127
0128
0129
0130 void compute(ResultType& res) const;
0131 };
0132
0133 template<typename MatrixType>
0134 MatrixPowerAtomic<MatrixType>::MatrixPowerAtomic(const MatrixType& T, RealScalar p) :
0135 m_A(T), m_p(p)
0136 {
0137 eigen_assert(T.rows() == T.cols());
0138 eigen_assert(p > -1 && p < 1);
0139 }
0140
0141 template<typename MatrixType>
0142 void MatrixPowerAtomic<MatrixType>::compute(ResultType& res) const
0143 {
0144 using std::pow;
0145 switch (m_A.rows()) {
0146 case 0:
0147 break;
0148 case 1:
0149 res(0,0) = pow(m_A(0,0), m_p);
0150 break;
0151 case 2:
0152 compute2x2(res, m_p);
0153 break;
0154 default:
0155 computeBig(res);
0156 }
0157 }
0158
0159 template<typename MatrixType>
0160 void MatrixPowerAtomic<MatrixType>::computePade(int degree, const MatrixType& IminusT, ResultType& res) const
0161 {
0162 int i = 2*degree;
0163 res = (m_p-RealScalar(degree)) / RealScalar(2*i-2) * IminusT;
0164
0165 for (--i; i; --i) {
0166 res = (MatrixType::Identity(IminusT.rows(), IminusT.cols()) + res).template triangularView<Upper>()
0167 .solve((i==1 ? -m_p : i&1 ? (-m_p-RealScalar(i/2))/RealScalar(2*i) : (m_p-RealScalar(i/2))/RealScalar(2*i-2)) * IminusT).eval();
0168 }
0169 res += MatrixType::Identity(IminusT.rows(), IminusT.cols());
0170 }
0171
0172
0173 template<typename MatrixType>
0174 void MatrixPowerAtomic<MatrixType>::compute2x2(ResultType& res, RealScalar p) const
0175 {
0176 using std::abs;
0177 using std::pow;
0178 res.coeffRef(0,0) = pow(m_A.coeff(0,0), p);
0179
0180 for (Index i=1; i < m_A.cols(); ++i) {
0181 res.coeffRef(i,i) = pow(m_A.coeff(i,i), p);
0182 if (m_A.coeff(i-1,i-1) == m_A.coeff(i,i))
0183 res.coeffRef(i-1,i) = p * pow(m_A.coeff(i,i), p-1);
0184 else if (2*abs(m_A.coeff(i-1,i-1)) < abs(m_A.coeff(i,i)) || 2*abs(m_A.coeff(i,i)) < abs(m_A.coeff(i-1,i-1)))
0185 res.coeffRef(i-1,i) = (res.coeff(i,i)-res.coeff(i-1,i-1)) / (m_A.coeff(i,i)-m_A.coeff(i-1,i-1));
0186 else
0187 res.coeffRef(i-1,i) = computeSuperDiag(m_A.coeff(i,i), m_A.coeff(i-1,i-1), p);
0188 res.coeffRef(i-1,i) *= m_A.coeff(i-1,i);
0189 }
0190 }
0191
0192 template<typename MatrixType>
0193 void MatrixPowerAtomic<MatrixType>::computeBig(ResultType& res) const
0194 {
0195 using std::ldexp;
0196 const int digits = std::numeric_limits<RealScalar>::digits;
0197 const RealScalar maxNormForPade = RealScalar(
0198 digits <= 24? 4.3386528e-1L
0199 : digits <= 53? 2.789358995219730e-1L
0200 : digits <= 64? 2.4471944416607995472e-1L
0201 : digits <= 106? 1.1016843812851143391275867258512e-1L
0202 : 9.134603732914548552537150753385375e-2L);
0203 MatrixType IminusT, sqrtT, T = m_A.template triangularView<Upper>();
0204 RealScalar normIminusT;
0205 int degree, degree2, numberOfSquareRoots = 0;
0206 bool hasExtraSquareRoot = false;
0207
0208 for (Index i=0; i < m_A.cols(); ++i)
0209 eigen_assert(m_A(i,i) != RealScalar(0));
0210
0211 while (true) {
0212 IminusT = MatrixType::Identity(m_A.rows(), m_A.cols()) - T;
0213 normIminusT = IminusT.cwiseAbs().colwise().sum().maxCoeff();
0214 if (normIminusT < maxNormForPade) {
0215 degree = getPadeDegree(normIminusT);
0216 degree2 = getPadeDegree(normIminusT/2);
0217 if (degree - degree2 <= 1 || hasExtraSquareRoot)
0218 break;
0219 hasExtraSquareRoot = true;
0220 }
0221 matrix_sqrt_triangular(T, sqrtT);
0222 T = sqrtT.template triangularView<Upper>();
0223 ++numberOfSquareRoots;
0224 }
0225 computePade(degree, IminusT, res);
0226
0227 for (; numberOfSquareRoots; --numberOfSquareRoots) {
0228 compute2x2(res, ldexp(m_p, -numberOfSquareRoots));
0229 res = res.template triangularView<Upper>() * res;
0230 }
0231 compute2x2(res, m_p);
0232 }
0233
0234 template<typename MatrixType>
0235 inline int MatrixPowerAtomic<MatrixType>::getPadeDegree(float normIminusT)
0236 {
0237 const float maxNormForPade[] = { 2.8064004e-1f , 4.3386528e-1f };
0238 int degree = 3;
0239 for (; degree <= 4; ++degree)
0240 if (normIminusT <= maxNormForPade[degree - 3])
0241 break;
0242 return degree;
0243 }
0244
0245 template<typename MatrixType>
0246 inline int MatrixPowerAtomic<MatrixType>::getPadeDegree(double normIminusT)
0247 {
0248 const double maxNormForPade[] = { 1.884160592658218e-2 , 6.038881904059573e-2, 1.239917516308172e-1,
0249 1.999045567181744e-1, 2.789358995219730e-1 };
0250 int degree = 3;
0251 for (; degree <= 7; ++degree)
0252 if (normIminusT <= maxNormForPade[degree - 3])
0253 break;
0254 return degree;
0255 }
0256
0257 template<typename MatrixType>
0258 inline int MatrixPowerAtomic<MatrixType>::getPadeDegree(long double normIminusT)
0259 {
0260 #if LDBL_MANT_DIG == 53
0261 const int maxPadeDegree = 7;
0262 const double maxNormForPade[] = { 1.884160592658218e-2L , 6.038881904059573e-2L, 1.239917516308172e-1L,
0263 1.999045567181744e-1L, 2.789358995219730e-1L };
0264 #elif LDBL_MANT_DIG <= 64
0265 const int maxPadeDegree = 8;
0266 const long double maxNormForPade[] = { 6.3854693117491799460e-3L , 2.6394893435456973676e-2L,
0267 6.4216043030404063729e-2L, 1.1701165502926694307e-1L, 1.7904284231268670284e-1L, 2.4471944416607995472e-1L };
0268 #elif LDBL_MANT_DIG <= 106
0269 const int maxPadeDegree = 10;
0270 const double maxNormForPade[] = { 1.0007161601787493236741409687186e-4L ,
0271 1.0007161601787493236741409687186e-3L, 4.7069769360887572939882574746264e-3L, 1.3220386624169159689406653101695e-2L,
0272 2.8063482381631737920612944054906e-2L, 4.9625993951953473052385361085058e-2L, 7.7367040706027886224557538328171e-2L,
0273 1.1016843812851143391275867258512e-1L };
0274 #else
0275 const int maxPadeDegree = 10;
0276 const double maxNormForPade[] = { 5.524506147036624377378713555116378e-5L ,
0277 6.640600568157479679823602193345995e-4L, 3.227716520106894279249709728084626e-3L,
0278 9.619593944683432960546978734646284e-3L, 2.134595382433742403911124458161147e-2L,
0279 3.908166513900489428442993794761185e-2L, 6.266780814639442865832535460550138e-2L,
0280 9.134603732914548552537150753385375e-2L };
0281 #endif
0282 int degree = 3;
0283 for (; degree <= maxPadeDegree; ++degree)
0284 if (normIminusT <= maxNormForPade[degree - 3])
0285 break;
0286 return degree;
0287 }
0288
0289 template<typename MatrixType>
0290 inline typename MatrixPowerAtomic<MatrixType>::ComplexScalar
0291 MatrixPowerAtomic<MatrixType>::computeSuperDiag(const ComplexScalar& curr, const ComplexScalar& prev, RealScalar p)
0292 {
0293 using std::ceil;
0294 using std::exp;
0295 using std::log;
0296 using std::sinh;
0297
0298 ComplexScalar logCurr = log(curr);
0299 ComplexScalar logPrev = log(prev);
0300 RealScalar unwindingNumber = ceil((numext::imag(logCurr - logPrev) - RealScalar(EIGEN_PI)) / RealScalar(2*EIGEN_PI));
0301 ComplexScalar w = numext::log1p((curr-prev)/prev)/RealScalar(2) + ComplexScalar(0, RealScalar(EIGEN_PI)*unwindingNumber);
0302 return RealScalar(2) * exp(RealScalar(0.5) * p * (logCurr + logPrev)) * sinh(p * w) / (curr - prev);
0303 }
0304
0305 template<typename MatrixType>
0306 inline typename MatrixPowerAtomic<MatrixType>::RealScalar
0307 MatrixPowerAtomic<MatrixType>::computeSuperDiag(RealScalar curr, RealScalar prev, RealScalar p)
0308 {
0309 using std::exp;
0310 using std::log;
0311 using std::sinh;
0312
0313 RealScalar w = numext::log1p((curr-prev)/prev)/RealScalar(2);
0314 return 2 * exp(p * (log(curr) + log(prev)) / 2) * sinh(p * w) / (curr - prev);
0315 }
0316
0317
0318
0319
0320
0321
0322
0323
0324
0325
0326
0327
0328
0329
0330
0331
0332
0333
0334
0335
0336 template<typename MatrixType>
0337 class MatrixPower : internal::noncopyable
0338 {
0339 private:
0340 typedef typename MatrixType::Scalar Scalar;
0341 typedef typename MatrixType::RealScalar RealScalar;
0342
0343 public:
0344
0345
0346
0347
0348
0349
0350
0351
0352 explicit MatrixPower(const MatrixType& A) :
0353 m_A(A),
0354 m_conditionNumber(0),
0355 m_rank(A.cols()),
0356 m_nulls(0)
0357 { eigen_assert(A.rows() == A.cols()); }
0358
0359
0360
0361
0362
0363
0364
0365
0366 const MatrixPowerParenthesesReturnValue<MatrixType> operator()(RealScalar p)
0367 { return MatrixPowerParenthesesReturnValue<MatrixType>(*this, p); }
0368
0369
0370
0371
0372
0373
0374
0375
0376 template<typename ResultType>
0377 void compute(ResultType& res, RealScalar p);
0378
0379 Index rows() const { return m_A.rows(); }
0380 Index cols() const { return m_A.cols(); }
0381
0382 private:
0383 typedef std::complex<RealScalar> ComplexScalar;
0384 typedef Matrix<ComplexScalar, Dynamic, Dynamic, 0,
0385 MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime> ComplexMatrix;
0386
0387
0388 typename MatrixType::Nested m_A;
0389
0390
0391 MatrixType m_tmp;
0392
0393
0394 ComplexMatrix m_T, m_U;
0395
0396
0397 ComplexMatrix m_fT;
0398
0399
0400
0401
0402
0403
0404
0405 RealScalar m_conditionNumber;
0406
0407
0408 Index m_rank;
0409
0410
0411 Index m_nulls;
0412
0413
0414
0415
0416
0417
0418
0419
0420
0421
0422 void split(RealScalar& p, RealScalar& intpart);
0423
0424
0425 void initialize();
0426
0427 template<typename ResultType>
0428 void computeIntPower(ResultType& res, RealScalar p);
0429
0430 template<typename ResultType>
0431 void computeFracPower(ResultType& res, RealScalar p);
0432
0433 template<int Rows, int Cols, int Options, int MaxRows, int MaxCols>
0434 static void revertSchur(
0435 Matrix<ComplexScalar, Rows, Cols, Options, MaxRows, MaxCols>& res,
0436 const ComplexMatrix& T,
0437 const ComplexMatrix& U);
0438
0439 template<int Rows, int Cols, int Options, int MaxRows, int MaxCols>
0440 static void revertSchur(
0441 Matrix<RealScalar, Rows, Cols, Options, MaxRows, MaxCols>& res,
0442 const ComplexMatrix& T,
0443 const ComplexMatrix& U);
0444 };
0445
0446 template<typename MatrixType>
0447 template<typename ResultType>
0448 void MatrixPower<MatrixType>::compute(ResultType& res, RealScalar p)
0449 {
0450 using std::pow;
0451 switch (cols()) {
0452 case 0:
0453 break;
0454 case 1:
0455 res(0,0) = pow(m_A.coeff(0,0), p);
0456 break;
0457 default:
0458 RealScalar intpart;
0459 split(p, intpart);
0460
0461 res = MatrixType::Identity(rows(), cols());
0462 computeIntPower(res, intpart);
0463 if (p) computeFracPower(res, p);
0464 }
0465 }
0466
0467 template<typename MatrixType>
0468 void MatrixPower<MatrixType>::split(RealScalar& p, RealScalar& intpart)
0469 {
0470 using std::floor;
0471 using std::pow;
0472
0473 intpart = floor(p);
0474 p -= intpart;
0475
0476
0477
0478 if (!m_conditionNumber && p)
0479 initialize();
0480
0481
0482 if (p > RealScalar(0.5) && p > (1-p) * pow(m_conditionNumber, p)) {
0483 --p;
0484 ++intpart;
0485 }
0486 }
0487
0488 template<typename MatrixType>
0489 void MatrixPower<MatrixType>::initialize()
0490 {
0491 const ComplexSchur<MatrixType> schurOfA(m_A);
0492 JacobiRotation<ComplexScalar> rot;
0493 ComplexScalar eigenvalue;
0494
0495 m_fT.resizeLike(m_A);
0496 m_T = schurOfA.matrixT();
0497 m_U = schurOfA.matrixU();
0498 m_conditionNumber = m_T.diagonal().array().abs().maxCoeff() / m_T.diagonal().array().abs().minCoeff();
0499
0500
0501 for (Index i = cols()-1; i>=0; --i) {
0502 if (m_rank <= 2)
0503 return;
0504 if (m_T.coeff(i,i) == RealScalar(0)) {
0505 for (Index j=i+1; j < m_rank; ++j) {
0506 eigenvalue = m_T.coeff(j,j);
0507 rot.makeGivens(m_T.coeff(j-1,j), eigenvalue);
0508 m_T.applyOnTheRight(j-1, j, rot);
0509 m_T.applyOnTheLeft(j-1, j, rot.adjoint());
0510 m_T.coeffRef(j-1,j-1) = eigenvalue;
0511 m_T.coeffRef(j,j) = RealScalar(0);
0512 m_U.applyOnTheRight(j-1, j, rot);
0513 }
0514 --m_rank;
0515 }
0516 }
0517
0518 m_nulls = rows() - m_rank;
0519 if (m_nulls) {
0520 eigen_assert(m_T.bottomRightCorner(m_nulls, m_nulls).isZero()
0521 && "Base of matrix power should be invertible or with a semisimple zero eigenvalue.");
0522 m_fT.bottomRows(m_nulls).fill(RealScalar(0));
0523 }
0524 }
0525
0526 template<typename MatrixType>
0527 template<typename ResultType>
0528 void MatrixPower<MatrixType>::computeIntPower(ResultType& res, RealScalar p)
0529 {
0530 using std::abs;
0531 using std::fmod;
0532 RealScalar pp = abs(p);
0533
0534 if (p<0)
0535 m_tmp = m_A.inverse();
0536 else
0537 m_tmp = m_A;
0538
0539 while (true) {
0540 if (fmod(pp, 2) >= 1)
0541 res = m_tmp * res;
0542 pp /= 2;
0543 if (pp < 1)
0544 break;
0545 m_tmp *= m_tmp;
0546 }
0547 }
0548
0549 template<typename MatrixType>
0550 template<typename ResultType>
0551 void MatrixPower<MatrixType>::computeFracPower(ResultType& res, RealScalar p)
0552 {
0553 Block<ComplexMatrix,Dynamic,Dynamic> blockTp(m_fT, 0, 0, m_rank, m_rank);
0554 eigen_assert(m_conditionNumber);
0555 eigen_assert(m_rank + m_nulls == rows());
0556
0557 MatrixPowerAtomic<ComplexMatrix>(m_T.topLeftCorner(m_rank, m_rank), p).compute(blockTp);
0558 if (m_nulls) {
0559 m_fT.topRightCorner(m_rank, m_nulls) = m_T.topLeftCorner(m_rank, m_rank).template triangularView<Upper>()
0560 .solve(blockTp * m_T.topRightCorner(m_rank, m_nulls));
0561 }
0562 revertSchur(m_tmp, m_fT, m_U);
0563 res = m_tmp * res;
0564 }
0565
0566 template<typename MatrixType>
0567 template<int Rows, int Cols, int Options, int MaxRows, int MaxCols>
0568 inline void MatrixPower<MatrixType>::revertSchur(
0569 Matrix<ComplexScalar, Rows, Cols, Options, MaxRows, MaxCols>& res,
0570 const ComplexMatrix& T,
0571 const ComplexMatrix& U)
0572 { res.noalias() = U * (T.template triangularView<Upper>() * U.adjoint()); }
0573
0574 template<typename MatrixType>
0575 template<int Rows, int Cols, int Options, int MaxRows, int MaxCols>
0576 inline void MatrixPower<MatrixType>::revertSchur(
0577 Matrix<RealScalar, Rows, Cols, Options, MaxRows, MaxCols>& res,
0578 const ComplexMatrix& T,
0579 const ComplexMatrix& U)
0580 { res.noalias() = (U * (T.template triangularView<Upper>() * U.adjoint())).real(); }
0581
0582
0583
0584
0585
0586
0587
0588
0589
0590
0591
0592
0593
0594
0595 template<typename Derived>
0596 class MatrixPowerReturnValue : public ReturnByValue< MatrixPowerReturnValue<Derived> >
0597 {
0598 public:
0599 typedef typename Derived::PlainObject PlainObject;
0600 typedef typename Derived::RealScalar RealScalar;
0601
0602
0603
0604
0605
0606
0607
0608 MatrixPowerReturnValue(const Derived& A, RealScalar p) : m_A(A), m_p(p)
0609 { }
0610
0611
0612
0613
0614
0615
0616
0617 template<typename ResultType>
0618 inline void evalTo(ResultType& result) const
0619 { MatrixPower<PlainObject>(m_A.eval()).compute(result, m_p); }
0620
0621 Index rows() const { return m_A.rows(); }
0622 Index cols() const { return m_A.cols(); }
0623
0624 private:
0625 const Derived& m_A;
0626 const RealScalar m_p;
0627 };
0628
0629
0630
0631
0632
0633
0634
0635
0636
0637
0638
0639
0640
0641
0642 template<typename Derived>
0643 class MatrixComplexPowerReturnValue : public ReturnByValue< MatrixComplexPowerReturnValue<Derived> >
0644 {
0645 public:
0646 typedef typename Derived::PlainObject PlainObject;
0647 typedef typename std::complex<typename Derived::RealScalar> ComplexScalar;
0648
0649
0650
0651
0652
0653
0654
0655 MatrixComplexPowerReturnValue(const Derived& A, const ComplexScalar& p) : m_A(A), m_p(p)
0656 { }
0657
0658
0659
0660
0661
0662
0663
0664
0665
0666
0667 template<typename ResultType>
0668 inline void evalTo(ResultType& result) const
0669 { result = (m_p * m_A.log()).exp(); }
0670
0671 Index rows() const { return m_A.rows(); }
0672 Index cols() const { return m_A.cols(); }
0673
0674 private:
0675 const Derived& m_A;
0676 const ComplexScalar m_p;
0677 };
0678
0679 namespace internal {
0680
0681 template<typename MatrixPowerType>
0682 struct traits< MatrixPowerParenthesesReturnValue<MatrixPowerType> >
0683 { typedef typename MatrixPowerType::PlainObject ReturnType; };
0684
0685 template<typename Derived>
0686 struct traits< MatrixPowerReturnValue<Derived> >
0687 { typedef typename Derived::PlainObject ReturnType; };
0688
0689 template<typename Derived>
0690 struct traits< MatrixComplexPowerReturnValue<Derived> >
0691 { typedef typename Derived::PlainObject ReturnType; };
0692
0693 }
0694
0695 template<typename Derived>
0696 const MatrixPowerReturnValue<Derived> MatrixBase<Derived>::pow(const RealScalar& p) const
0697 { return MatrixPowerReturnValue<Derived>(derived(), p); }
0698
0699 template<typename Derived>
0700 const MatrixComplexPowerReturnValue<Derived> MatrixBase<Derived>::pow(const std::complex<RealScalar>& p) const
0701 { return MatrixComplexPowerReturnValue<Derived>(derived(), p); }
0702
0703 }
0704
0705 #endif