10#ifndef EIGEN_MATRIX_FUNCTION_H
11#define EIGEN_MATRIX_FUNCTION_H
31template <
typename MatrixType>
34 typedef typename MatrixType::Scalar
Scalar;
46 MatrixType
compute(
const MatrixType& A);
52template <
typename MatrixType>
55 Index rows = A.rows();
56 const MatrixType N = MatrixType::Identity(rows, rows) - A;
57 VectorType
e = VectorType::Ones(rows);
58 N.template triangularView<Upper>().solveInPlace(
e);
59 return e.cwiseAbs().maxCoeff();
62template <
typename MatrixType>
65 typedef typename NumTraits<Scalar>::Real RealScalar;
66 Index rows = A.rows();
68 MatrixType Ashifted = A - avgEival * MatrixType::Identity(rows, rows);
70 MatrixType
F = m_f(avgEival, 0) * MatrixType::Identity(rows, rows);
71 MatrixType P = Ashifted;
73 for (Index s = 1; double(s) < 1.1 * double(rows) + 10.0; s++) {
74 Fincr = m_f(avgEival,
static_cast<int>(s)) * P;
76 P =
Scalar(RealScalar(1) / RealScalar(s + 1)) * P * Ashifted;
79 const RealScalar F_norm =
F.cwiseAbs().rowwise().sum().maxCoeff();
80 const RealScalar Fincr_norm = Fincr.cwiseAbs().rowwise().sum().maxCoeff();
81 if (Fincr_norm < NumTraits<Scalar>::epsilon() * F_norm) {
83 RealScalar rfactorial = 1;
84 for (Index r = 0; r < rows; r++) {
86 for (Index i = 0; i < rows; i++)
87 mx = (
std::max)(mx,
std::abs(m_f(Ashifted(i, i) + avgEival,
static_cast<int>(s + r))));
88 if (r != 0) rfactorial *= RealScalar(r);
89 delta = (
std::max)(delta, mx / rfactorial);
91 const RealScalar P_norm = P.cwiseAbs().rowwise().sum().maxCoeff();
92 if (mu * delta * P_norm < NumTraits<Scalar>::epsilon() * F_norm)
104template <
typename Index,
typename ListOfClusters>
106 typename std::list<Index>::iterator j;
107 for (
typename ListOfClusters::iterator i = clusters.begin(); i != clusters.end(); ++i) {
108 j =
std::find(i->begin(), i->end(), key);
109 if (j != i->end())
return i;
111 return clusters.end();
125template <
typename EivalsType,
typename Cluster>
127 typedef typename EivalsType::RealScalar RealScalar;
128 for (Index i = 0; i < eivals.rows(); ++i) {
131 if (qi == clusters.end()) {
134 clusters.push_back(l);
140 for (Index j = i + 1; j < eivals.rows(); ++j) {
142 std::find(qi->begin(), qi->end(), j) == qi->end()) {
144 if (qj == clusters.end()) {
147 qi->insert(qi->end(), qj->begin(), qj->end());
156template <
typename ListOfClusters,
typename Index>
158 const Index numClusters =
static_cast<Index
>(clusters.size());
159 clusterSize.setZero(numClusters);
160 Index clusterIndex = 0;
161 for (
typename ListOfClusters::const_iterator cluster = clusters.begin(); cluster != clusters.end(); ++cluster) {
162 clusterSize[clusterIndex] = cluster->size();
168template <
typename VectorType>
170 blockStart.resize(clusterSize.rows());
172 for (Index i = 1; i < clusterSize.rows(); i++) {
173 blockStart(i) = blockStart(i - 1) + clusterSize(i - 1);
178template <
typename EivalsType,
typename ListOfClusters,
typename VectorType>
180 eivalToCluster.resize(eivals.rows());
181 Index clusterIndex = 0;
182 for (
typename ListOfClusters::const_iterator cluster = clusters.begin(); cluster != clusters.end(); ++cluster) {
183 for (Index i = 0; i < eivals.rows(); ++i) {
184 if (
std::find(cluster->begin(), cluster->end(), i) != cluster->end()) {
185 eivalToCluster[i] = clusterIndex;
193template <
typename DynVectorType,
typename VectorType>
195 VectorType& permutation) {
196 DynVectorType indexNextEntry = blockStart;
197 permutation.resize(eivalToCluster.rows());
198 for (Index i = 0; i < eivalToCluster.rows(); i++) {
199 Index cluster = eivalToCluster[i];
200 permutation[i] = indexNextEntry[cluster];
201 ++indexNextEntry[cluster];
206template <
typename VectorType,
typename MatrixType>
208 for (Index i = 0; i < permutation.rows() - 1; i++) {
210 for (j = i; j < permutation.rows(); j++) {
211 if (permutation(j) == i)
break;
213 eigen_assert(permutation(j) == i);
214 for (Index k = j - 1; k >= i; k--) {
215 JacobiRotation<typename MatrixType::Scalar> rotation;
216 rotation.makeGivens(T(k, k + 1), T(k + 1, k + 1) - T(k, k));
217 T.applyOnTheLeft(k, k + 1, rotation.adjoint());
218 T.applyOnTheRight(k, k + 1, rotation);
219 U.applyOnTheRight(k, k + 1, rotation);
220 std::swap(permutation.coeffRef(k), permutation.coeffRef(k + 1));
231template <
typename MatrixType,
typename AtomicType,
typename VectorType>
233 const VectorType& clusterSize, MatrixType& fT) {
234 fT.setZero(T.rows(), T.cols());
235 for (Index i = 0; i < clusterSize.rows(); ++i) {
236 fT.block(blockStart(i), blockStart(i), clusterSize(i), clusterSize(i)) =
237 atomic.compute(T.block(blockStart(i), blockStart(i), clusterSize(i), clusterSize(i)));
263template <
typename MatrixType>
265 eigen_assert(A.rows() == A.cols());
266 eigen_assert(A.isUpperTriangular());
267 eigen_assert(B.rows() == B.cols());
268 eigen_assert(B.isUpperTriangular());
269 eigen_assert(C.rows() == A.rows());
270 eigen_assert(C.cols() == B.rows());
272 typedef typename MatrixType::Scalar Scalar;
278 for (Index i = m - 1; i >= 0; --i) {
279 for (Index j = 0; j < n; ++j) {
285 Matrix<Scalar, 1, 1> AXmatrix = A.row(i).tail(m - 1 - i) * X.col(j).tail(m - 1 - i);
294 Matrix<Scalar, 1, 1> XBmatrix = X.row(i).head(j) * B.col(j).head(j);
298 X(i, j) = (C(i, j) - AX - XB) / (A(i, i) + B(j, j));
310template <
typename MatrixType,
typename VectorType>
312 const VectorType& clusterSize, MatrixType& fT) {
313 typedef internal::traits<MatrixType> Traits;
314 typedef typename MatrixType::Scalar Scalar;
315 static const int Options = MatrixType::Options;
316 typedef Matrix<Scalar, Dynamic, Dynamic, Options, Traits::RowsAtCompileTime, Traits::ColsAtCompileTime> DynMatrixType;
318 for (Index k = 1; k < clusterSize.rows(); k++) {
319 for (Index i = 0; i < clusterSize.rows() - k; i++) {
321 DynMatrixType A = T.block(blockStart(i), blockStart(i), clusterSize(i), clusterSize(i));
322 DynMatrixType B = -T.block(blockStart(i + k), blockStart(i + k), clusterSize(i + k), clusterSize(i + k));
323 DynMatrixType C = fT.block(blockStart(i), blockStart(i), clusterSize(i), clusterSize(i)) *
324 T.block(blockStart(i), blockStart(i + k), clusterSize(i), clusterSize(i + k));
325 C -= T.block(blockStart(i), blockStart(i + k), clusterSize(i), clusterSize(i + k)) *
326 fT.block(blockStart(i + k), blockStart(i + k), clusterSize(i + k), clusterSize(i + k));
327 for (Index m = i + 1; m < i + k; m++) {
328 C += fT.block(blockStart(i), blockStart(m), clusterSize(i), clusterSize(m)) *
329 T.block(blockStart(m), blockStart(i + k), clusterSize(m), clusterSize(i + k));
330 C -= T.block(blockStart(i), blockStart(m), clusterSize(i), clusterSize(m)) *
331 fT.block(blockStart(m), blockStart(i + k), clusterSize(m), clusterSize(i + k));
333 fT.block(blockStart(i), blockStart(i + k), clusterSize(i), clusterSize(i + k)) =
354template <typename MatrixType, int IsComplex = NumTraits<typename internal::traits<MatrixType>::Scalar>::IsComplex>
366 template <
typename AtomicType,
typename ResultType>
367 static void run(
const MatrixType& A, AtomicType& atomic, ResultType& result);
376template <
typename MatrixType>
378 template <
typename MatA,
typename AtomicType,
typename ResultType>
379 static void run(
const MatA& A, AtomicType& atomic, ResultType& result) {
380 typedef internal::traits<MatrixType> Traits;
381 typedef typename Traits::Scalar Scalar;
382 static const int Rows = Traits::RowsAtCompileTime, Cols = Traits::ColsAtCompileTime;
383 static const int MaxRows = Traits::MaxRowsAtCompileTime, MaxCols = Traits::MaxColsAtCompileTime;
385 typedef std::complex<Scalar> ComplexScalar;
386 typedef Matrix<ComplexScalar, Rows, Cols, 0, MaxRows, MaxCols> ComplexMatrix;
388 ComplexMatrix CA = A.template cast<ComplexScalar>();
389 ComplexMatrix Cresult;
391 result = Cresult.real();
398template <
typename MatrixType>
400 template <
typename MatA,
typename AtomicType,
typename ResultType>
401 static void run(
const MatA& A, AtomicType& atomic, ResultType& result) {
402 typedef internal::traits<MatrixType> Traits;
405 const ComplexSchur<MatrixType> schurOfA(A);
406 eigen_assert(schurOfA.info() == Success);
407 MatrixType T = schurOfA.matrixT();
408 MatrixType U = schurOfA.matrixU();
411 std::list<std::list<Index> > clusters;
415 Matrix<Index, Dynamic, 1> clusterSize;
419 Matrix<Index, Dynamic, 1> blockStart;
423 Matrix<Index, Dynamic, 1> eivalToCluster;
427 Matrix<Index, Traits::RowsAtCompileTime, 1> permutation;
437 result = U * (fT.template triangularView<Upper>() * U.adjoint());
453template <
typename Derived>
474 template <
typename ResultType>
475 inline void evalTo(ResultType& result)
const {
477 typedef internal::remove_all_t<NestedEvalType> NestedEvalTypeClean;
478 typedef internal::traits<NestedEvalTypeClean> Traits;
479 typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
480 typedef Matrix<ComplexScalar, Dynamic, Dynamic, 0, Traits::RowsAtCompileTime, Traits::ColsAtCompileTime>
484 AtomicType atomic(m_f);
489 Index
rows()
const {
return m_A.rows(); }
490 Index
cols()
const {
return m_A.cols(); }
498template <
typename Derived>
506template <
typename Derived>
508 typename internal::stem_function<
typename internal::traits<Derived>::Scalar>
::type f)
const {
509 eigen_assert(rows() == cols());
513template <
typename Derived>
515 eigen_assert(rows() == cols());
516 typedef typename internal::stem_function<Scalar>::ComplexScalar ComplexScalar;
517 return MatrixFunctionReturnValue<Derived>(derived(), internal::stem_function_sin<ComplexScalar>);
520template <
typename Derived>
522 eigen_assert(rows() == cols());
523 typedef typename internal::stem_function<Scalar>::ComplexScalar ComplexScalar;
524 return MatrixFunctionReturnValue<Derived>(derived(), internal::stem_function_cos<ComplexScalar>);
527template <
typename Derived>
529 eigen_assert(rows() == cols());
530 typedef typename internal::stem_function<Scalar>::ComplexScalar ComplexScalar;
531 return MatrixFunctionReturnValue<Derived>(derived(), internal::stem_function_sinh<ComplexScalar>);
534template <
typename Derived>
536 eigen_assert(rows() == cols());
537 typedef typename internal::stem_function<Scalar>::ComplexScalar ComplexScalar;
538 return MatrixFunctionReturnValue<Derived>(derived(), internal::stem_function_cosh<ComplexScalar>);
Proxy for the matrix function of some matrix (expression).
Definition: MatrixFunction.h:454
internal::ref_selector< Derived >::type DerivedNested
Definition: MatrixFunction.h:460
Index rows() const
Definition: MatrixFunction.h:489
internal::stem_function< Scalar >::type StemFunction
Definition: MatrixFunction.h:457
void evalTo(ResultType &result) const
Compute the matrix function.
Definition: MatrixFunction.h:475
Index cols() const
Definition: MatrixFunction.h:490
Derived::Scalar Scalar
Definition: MatrixFunction.h:456
MatrixFunctionReturnValue(const Derived &A, StemFunction f)
Constructor.
Definition: MatrixFunction.h:468
Helper class for computing matrix functions of atomic matrices.
Definition: MatrixFunction.h:32
MatrixFunctionAtomic(StemFunction f)
Constructor.
Definition: MatrixFunction.h:40
MatrixType::Scalar Scalar
Definition: MatrixFunction.h:34
MatrixType compute(const MatrixType &A)
Compute matrix function of atomic matrix.
Definition: MatrixFunction.h:63
stem_function< Scalar >::type StemFunction
Definition: MatrixFunction.h:35
UnitType abs(const UnitType x) noexcept
Compute absolute value.
Definition: math.h:721
dimensionless::scalar_t sinh(const AngleUnit angle) noexcept
Compute hyperbolic sine.
Definition: math.h:226
dimensionless::scalar_t cosh(const AngleUnit angle) noexcept
Compute hyperbolic cosine.
Definition: math.h:206
dimensionless::scalar_t cos(const AngleUnit angle) noexcept
Compute cosine.
Definition: math.h:61
dimensionless::scalar_t sin(const AngleUnit angle) noexcept
Compute sine.
Definition: math.h:81
void matrix_function_compute_permutation(const DynVectorType &blockStart, const DynVectorType &eivalToCluster, VectorType &permutation)
Compute permutation which groups ei'vals in same cluster together.
Definition: MatrixFunction.h:194
void matrix_function_compute_cluster_size(const ListOfClusters &clusters, Matrix< Index, Dynamic, 1 > &clusterSize)
Compute size of each cluster given a partitioning.
Definition: MatrixFunction.h:157
void matrix_function_compute_block_start(const VectorType &clusterSize, VectorType &blockStart)
Compute start of each block using clusterSize.
Definition: MatrixFunction.h:169
void matrix_function_compute_block_atomic(const MatrixType &T, AtomicType &atomic, const VectorType &blockStart, const VectorType &clusterSize, MatrixType &fT)
Compute block diagonal part of matrix function.
Definition: MatrixFunction.h:232
void matrix_function_permute_schur(VectorType &permutation, MatrixType &U, MatrixType &T)
Permute Schur decomposition in U and T according to permutation.
Definition: MatrixFunction.h:207
static const float matrix_function_separation
Maximum distance allowed between eigenvalues to be considered "close".
Definition: MatrixFunction.h:23
void matrix_function_compute_above_diagonal(const MatrixType &T, const VectorType &blockStart, const VectorType &clusterSize, MatrixType &fT)
Compute part of matrix function above block diagonal.
Definition: MatrixFunction.h:311
void matrix_function_partition_eigenvalues(const EivalsType &eivals, std::list< Cluster > &clusters)
Partition eigenvalues in clusters of ei'vals close to each other.
Definition: MatrixFunction.h:126
MatrixType matrix_function_solve_triangular_sylvester(const MatrixType &A, const MatrixType &B, const MatrixType &C)
Solve a triangular Sylvester equation AX + XB = C.
Definition: MatrixFunction.h:264
void matrix_function_compute_map(const EivalsType &eivals, const ListOfClusters &clusters, VectorType &eivalToCluster)
Compute mapping of eigenvalue indices to cluster indices.
Definition: MatrixFunction.h:179
NumTraits< typenameMatrixType::Scalar >::Real matrix_function_compute_mu(const MatrixType &A)
Definition: MatrixFunction.h:53
ListOfClusters::iterator matrix_function_find_cluster(Index key, ListOfClusters &clusters)
Find cluster in clusters containing some value.
Definition: MatrixFunction.h:105
Definition: MatrixSquareRoot.h:16
FMT_CONSTEXPR auto find(Ptr first, Ptr last, T value, Ptr &out) -> bool
Definition: core.h:2120
type
Definition: core.h:556
WPI_BASIC_JSON_TPL_DECLARATION void swap(wpi::WPI_BASIC_JSON_TPL &j1, wpi::WPI_BASIC_JSON_TPL &j2) noexcept(//NOLINT(readability-inconsistent-declaration-parameter-name) is_nothrow_move_constructible< wpi::WPI_BASIC_JSON_TPL >::value &&//NOLINT(misc-redundant-expression) is_nothrow_move_assignable< wpi::WPI_BASIC_JSON_TPL >::value)
exchanges the values of two JSON objects
Definition: json.h:5219
static constexpr const unit_t< compound_unit< charge::coulomb, inverse< substance::mol > > > F(N_A *e)
Faraday constant.
static constexpr const charge::coulomb_t e(1.6021766208e-19)
elementary charge.
UnitTypeLhs() max(const UnitTypeLhs &lhs, const UnitTypeRhs &rhs)
Definition: base.h:3417
static void run(const MatA &A, AtomicType &atomic, ResultType &result)
Definition: MatrixFunction.h:379
static void run(const MatA &A, AtomicType &atomic, ResultType &result)
Definition: MatrixFunction.h:401
Class for computing matrix functions.
Definition: MatrixFunction.h:355
static void run(const MatrixType &A, AtomicType &atomic, ResultType &result)
Compute the matrix function.
Derived::PlainObject ReturnType
Definition: MatrixFunction.h:500