MatrixFunction.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) 2009-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_FUNCTION_H
11 #define EIGEN_MATRIX_FUNCTION_H
12 
13 #include "StemFunction.h"
14 
15 
16 namespace Eigen {
17 
18 namespace internal {
19 
21 static const float matrix_function_separation = 0.1f;
22 
29 template <typename MatrixType>
31 {
32  public:
33 
34  typedef typename MatrixType::Scalar Scalar;
36 
41 
47 
48  private:
50 };
51 
52 template <typename MatrixType>
54 {
56  Index rows = A.rows();
57  const MatrixType N = MatrixType::Identity(rows, rows) - A;
58  VectorType e = VectorType::Ones(rows);
59  N.template triangularView<Upper>().solveInPlace(e);
60  return e.cwiseAbs().maxCoeff();
61 }
62 
63 template <typename MatrixType>
65 {
66  // TODO: Use that A is upper triangular
67  typedef typename NumTraits<Scalar>::Real RealScalar;
68  Index rows = A.rows();
69  Scalar avgEival = A.trace() / Scalar(RealScalar(rows));
70  MatrixType Ashifted = A - avgEival * MatrixType::Identity(rows, rows);
72  MatrixType F = m_f(avgEival, 0) * MatrixType::Identity(rows, rows);
73  MatrixType P = Ashifted;
74  MatrixType Fincr;
75  for (Index s = 1; double(s) < 1.1 * double(rows) + 10.0; s++) { // upper limit is fairly arbitrary
76  Fincr = m_f(avgEival, static_cast<int>(s)) * P;
77  F += Fincr;
78  P = Scalar(RealScalar(1)/RealScalar(s + 1)) * P * Ashifted;
79 
80  // test whether Taylor series converged
81  const RealScalar F_norm = F.cwiseAbs().rowwise().sum().maxCoeff();
82  const RealScalar Fincr_norm = Fincr.cwiseAbs().rowwise().sum().maxCoeff();
83  if (Fincr_norm < NumTraits<Scalar>::epsilon() * F_norm) {
84  RealScalar delta = 0;
85  RealScalar rfactorial = 1;
86  for (Index r = 0; r < rows; r++) {
87  RealScalar mx = 0;
88  for (Index i = 0; i < rows; i++)
89  mx = (std::max)(mx, std::abs(m_f(Ashifted(i, i) + avgEival, static_cast<int>(s+r))));
90  if (r != 0)
91  rfactorial *= RealScalar(r);
92  delta = (std::max)(delta, mx / rfactorial);
93  }
94  const RealScalar P_norm = P.cwiseAbs().rowwise().sum().maxCoeff();
95  if (mu * delta * P_norm < NumTraits<Scalar>::epsilon() * F_norm) // series converged
96  break;
97  }
98  }
99  return F;
100 }
101 
107 template <typename Index, typename ListOfClusters>
108 typename ListOfClusters::iterator matrix_function_find_cluster(Index key, ListOfClusters& clusters)
109 {
110  typename std::list<Index>::iterator j;
111  for (typename ListOfClusters::iterator i = clusters.begin(); i != clusters.end(); ++i) {
112  j = std::find(i->begin(), i->end(), key);
113  if (j != i->end())
114  return i;
115  }
116  return clusters.end();
117 }
118 
130 template <typename EivalsType, typename Cluster>
131 void matrix_function_partition_eigenvalues(const EivalsType& eivals, std::list<Cluster>& clusters)
132 {
133  typedef typename EivalsType::RealScalar RealScalar;
134  for (Index i=0; i<eivals.rows(); ++i) {
135  // Find cluster containing i-th ei'val, adding a new cluster if necessary
136  typename std::list<Cluster>::iterator qi = matrix_function_find_cluster(i, clusters);
137  if (qi == clusters.end()) {
138  Cluster l;
139  l.push_back(i);
140  clusters.push_back(l);
141  qi = clusters.end();
142  --qi;
143  }
144 
145  // Look for other element to add to the set
146  for (Index j=i+1; j<eivals.rows(); ++j) {
148  && std::find(qi->begin(), qi->end(), j) == qi->end()) {
149  typename std::list<Cluster>::iterator qj = matrix_function_find_cluster(j, clusters);
150  if (qj == clusters.end()) {
151  qi->push_back(j);
152  } else {
153  qi->insert(qi->end(), qj->begin(), qj->end());
154  clusters.erase(qj);
155  }
156  }
157  }
158  }
159 }
160 
162 template <typename ListOfClusters, typename Index>
163 void matrix_function_compute_cluster_size(const ListOfClusters& clusters, Matrix<Index, Dynamic, 1>& clusterSize)
164 {
165  const Index numClusters = static_cast<Index>(clusters.size());
166  clusterSize.setZero(numClusters);
167  Index clusterIndex = 0;
168  for (typename ListOfClusters::const_iterator cluster = clusters.begin(); cluster != clusters.end(); ++cluster) {
169  clusterSize[clusterIndex] = cluster->size();
170  ++clusterIndex;
171  }
172 }
173 
175 template <typename VectorType>
176 void matrix_function_compute_block_start(const VectorType& clusterSize, VectorType& blockStart)
177 {
178  blockStart.resize(clusterSize.rows());
179  blockStart(0) = 0;
180  for (Index i = 1; i < clusterSize.rows(); i++) {
181  blockStart(i) = blockStart(i-1) + clusterSize(i-1);
182  }
183 }
184 
186 template <typename EivalsType, typename ListOfClusters, typename VectorType>
187 void matrix_function_compute_map(const EivalsType& eivals, const ListOfClusters& clusters, VectorType& eivalToCluster)
188 {
189  eivalToCluster.resize(eivals.rows());
190  Index clusterIndex = 0;
191  for (typename ListOfClusters::const_iterator cluster = clusters.begin(); cluster != clusters.end(); ++cluster) {
192  for (Index i = 0; i < eivals.rows(); ++i) {
193  if (std::find(cluster->begin(), cluster->end(), i) != cluster->end()) {
194  eivalToCluster[i] = clusterIndex;
195  }
196  }
197  ++clusterIndex;
198  }
199 }
200 
202 template <typename DynVectorType, typename VectorType>
203 void matrix_function_compute_permutation(const DynVectorType& blockStart, const DynVectorType& eivalToCluster, VectorType& permutation)
204 {
205  DynVectorType indexNextEntry = blockStart;
206  permutation.resize(eivalToCluster.rows());
207  for (Index i = 0; i < eivalToCluster.rows(); i++) {
208  Index cluster = eivalToCluster[i];
209  permutation[i] = indexNextEntry[cluster];
210  ++indexNextEntry[cluster];
211  }
212 }
213 
215 template <typename VectorType, typename MatrixType>
217 {
218  for (Index i = 0; i < permutation.rows() - 1; i++) {
219  Index j;
220  for (j = i; j < permutation.rows(); j++) {
221  if (permutation(j) == i) break;
222  }
223  eigen_assert(permutation(j) == i);
224  for (Index k = j-1; k >= i; k--) {
226  rotation.makeGivens(T(k, k+1), T(k+1, k+1) - T(k, k));
227  T.applyOnTheLeft(k, k+1, rotation.adjoint());
228  T.applyOnTheRight(k, k+1, rotation);
229  U.applyOnTheRight(k, k+1, rotation);
230  std::swap(permutation.coeffRef(k), permutation.coeffRef(k+1));
231  }
232  }
233 }
234 
241 template <typename MatrixType, typename AtomicType, typename VectorType>
242 void matrix_function_compute_block_atomic(const MatrixType& T, AtomicType& atomic, const VectorType& blockStart, const VectorType& clusterSize, MatrixType& fT)
243 {
244  fT.setZero(T.rows(), T.cols());
245  for (Index i = 0; i < clusterSize.rows(); ++i) {
246  fT.block(blockStart(i), blockStart(i), clusterSize(i), clusterSize(i))
247  = atomic.compute(T.block(blockStart(i), blockStart(i), clusterSize(i), clusterSize(i)));
248  }
249 }
250 
273 template <typename MatrixType>
275 {
276  eigen_assert(A.rows() == A.cols());
277  eigen_assert(A.isUpperTriangular());
278  eigen_assert(B.rows() == B.cols());
279  eigen_assert(B.isUpperTriangular());
280  eigen_assert(C.rows() == A.rows());
281  eigen_assert(C.cols() == B.rows());
282 
283  typedef typename MatrixType::Scalar Scalar;
284 
285  Index m = A.rows();
286  Index n = B.rows();
287  MatrixType X(m, n);
288 
289  for (Index i = m - 1; i >= 0; --i) {
290  for (Index j = 0; j < n; ++j) {
291 
292  // Compute AX = \sum_{k=i+1}^m A_{ik} X_{kj}
293  Scalar AX;
294  if (i == m - 1) {
295  AX = 0;
296  } else {
297  Matrix<Scalar,1,1> AXmatrix = A.row(i).tail(m-1-i) * X.col(j).tail(m-1-i);
298  AX = AXmatrix(0,0);
299  }
300 
301  // Compute XB = \sum_{k=1}^{j-1} X_{ik} B_{kj}
302  Scalar XB;
303  if (j == 0) {
304  XB = 0;
305  } else {
306  Matrix<Scalar,1,1> XBmatrix = X.row(i).head(j) * B.col(j).head(j);
307  XB = XBmatrix(0,0);
308  }
309 
310  X(i,j) = (C(i,j) - AX - XB) / (A(i,i) + B(j,j));
311  }
312  }
313  return X;
314 }
315 
322 template <typename MatrixType, typename VectorType>
323 void matrix_function_compute_above_diagonal(const MatrixType& T, const VectorType& blockStart, const VectorType& clusterSize, MatrixType& fT)
324 {
325  typedef internal::traits<MatrixType> Traits;
326  typedef typename MatrixType::Scalar Scalar;
327  static const int Options = MatrixType::Options;
329 
330  for (Index k = 1; k < clusterSize.rows(); k++) {
331  for (Index i = 0; i < clusterSize.rows() - k; i++) {
332  // compute (i, i+k) block
333  DynMatrixType A = T.block(blockStart(i), blockStart(i), clusterSize(i), clusterSize(i));
334  DynMatrixType B = -T.block(blockStart(i+k), blockStart(i+k), clusterSize(i+k), clusterSize(i+k));
335  DynMatrixType C = fT.block(blockStart(i), blockStart(i), clusterSize(i), clusterSize(i))
336  * T.block(blockStart(i), blockStart(i+k), clusterSize(i), clusterSize(i+k));
337  C -= T.block(blockStart(i), blockStart(i+k), clusterSize(i), clusterSize(i+k))
338  * fT.block(blockStart(i+k), blockStart(i+k), clusterSize(i+k), clusterSize(i+k));
339  for (Index m = i + 1; m < i + k; m++) {
340  C += fT.block(blockStart(i), blockStart(m), clusterSize(i), clusterSize(m))
341  * T.block(blockStart(m), blockStart(i+k), clusterSize(m), clusterSize(i+k));
342  C -= T.block(blockStart(i), blockStart(m), clusterSize(i), clusterSize(m))
343  * fT.block(blockStart(m), blockStart(i+k), clusterSize(m), clusterSize(i+k));
344  }
345  fT.block(blockStart(i), blockStart(i+k), clusterSize(i), clusterSize(i+k))
347  }
348  }
349 }
350 
368 {
379  template <typename AtomicType, typename ResultType>
380  static void run(const MatrixType& A, AtomicType& atomic, ResultType &result);
381 };
382 
389 template <typename MatrixType>
391 {
392  template <typename MatA, typename AtomicType, typename ResultType>
393  static void run(const MatA& A, AtomicType& atomic, ResultType &result)
394  {
395  typedef internal::traits<MatrixType> Traits;
396  typedef typename Traits::Scalar Scalar;
397  static const int Rows = Traits::RowsAtCompileTime, Cols = Traits::ColsAtCompileTime;
398  static const int MaxRows = Traits::MaxRowsAtCompileTime, MaxCols = Traits::MaxColsAtCompileTime;
399 
400  typedef std::complex<Scalar> ComplexScalar;
402 
403  ComplexMatrix CA = A.template cast<ComplexScalar>();
404  ComplexMatrix Cresult;
405  matrix_function_compute<ComplexMatrix>::run(CA, atomic, Cresult);
406  result = Cresult.real();
407  }
408 };
409 
413 template <typename MatrixType>
415 {
416  template <typename MatA, typename AtomicType, typename ResultType>
417  static void run(const MatA& A, AtomicType& atomic, ResultType &result)
418  {
419  typedef internal::traits<MatrixType> Traits;
420 
421  // compute Schur decomposition of A
423  eigen_assert(schurOfA.info()==Success);
424  MatrixType T = schurOfA.matrixT();
425  MatrixType U = schurOfA.matrixU();
426 
427  // partition eigenvalues into clusters of ei'vals "close" to each other
428  std::list<std::list<Index> > clusters;
429  matrix_function_partition_eigenvalues(T.diagonal(), clusters);
430 
431  // compute size of each cluster
432  Matrix<Index, Dynamic, 1> clusterSize;
433  matrix_function_compute_cluster_size(clusters, clusterSize);
434 
435  // blockStart[i] is row index at which block corresponding to i-th cluster starts
436  Matrix<Index, Dynamic, 1> blockStart;
437  matrix_function_compute_block_start(clusterSize, blockStart);
438 
439  // compute map so that eivalToCluster[i] = j means that i-th ei'val is in j-th cluster
440  Matrix<Index, Dynamic, 1> eivalToCluster;
441  matrix_function_compute_map(T.diagonal(), clusters, eivalToCluster);
442 
443  // compute permutation which groups ei'vals in same cluster together
445  matrix_function_compute_permutation(blockStart, eivalToCluster, permutation);
446 
447  // permute Schur decomposition
448  matrix_function_permute_schur(permutation, U, T);
449 
450  // compute result
451  MatrixType fT; // matrix function applied to T
452  matrix_function_compute_block_atomic(T, atomic, blockStart, clusterSize, fT);
453  matrix_function_compute_above_diagonal(T, blockStart, clusterSize, fT);
454  result = U * (fT.template triangularView<Upper>() * U.adjoint());
455  }
456 };
457 
458 } // end of namespace internal
459 
470 template<typename Derived> class MatrixFunctionReturnValue
471 : public ReturnByValue<MatrixFunctionReturnValue<Derived> >
472 {
473  public:
474  typedef typename Derived::Scalar Scalar;
476 
477  protected:
479 
480  public:
481 
487  MatrixFunctionReturnValue(const Derived& A, StemFunction f) : m_A(A), m_f(f) { }
488 
493  template <typename ResultType>
494  inline void evalTo(ResultType& result) const
495  {
496  typedef typename internal::nested_eval<Derived, 10>::type NestedEvalType;
497  typedef typename internal::remove_all<NestedEvalType>::type NestedEvalTypeClean;
499  typedef std::complex<typename NumTraits<Scalar>::Real> ComplexScalar;
501 
503  AtomicType atomic(m_f);
504 
506  }
507 
508  Index rows() const { return m_A.rows(); }
509  Index cols() const { return m_A.cols(); }
510 
511  private:
514 };
515 
516 namespace internal {
517 template<typename Derived>
519 {
520  typedef typename Derived::PlainObject ReturnType;
521 };
522 }
523 
524 
525 /********** MatrixBase methods **********/
526 
527 
528 template <typename Derived>
530 {
531  eigen_assert(rows() == cols());
532  return MatrixFunctionReturnValue<Derived>(derived(), f);
533 }
534 
535 template <typename Derived>
537 {
538  eigen_assert(rows() == cols());
539  typedef typename internal::stem_function<Scalar>::ComplexScalar ComplexScalar;
540  return MatrixFunctionReturnValue<Derived>(derived(), internal::stem_function_sin<ComplexScalar>);
541 }
542 
543 template <typename Derived>
544 const MatrixFunctionReturnValue<Derived> MatrixBase<Derived>::cos() const
545 {
546  eigen_assert(rows() == cols());
547  typedef typename internal::stem_function<Scalar>::ComplexScalar ComplexScalar;
548  return MatrixFunctionReturnValue<Derived>(derived(), internal::stem_function_cos<ComplexScalar>);
549 }
550 
551 template <typename Derived>
552 const MatrixFunctionReturnValue<Derived> MatrixBase<Derived>::sinh() const
553 {
554  eigen_assert(rows() == cols());
555  typedef typename internal::stem_function<Scalar>::ComplexScalar ComplexScalar;
556  return MatrixFunctionReturnValue<Derived>(derived(), internal::stem_function_sinh<ComplexScalar>);
557 }
558 
559 template <typename Derived>
560 const MatrixFunctionReturnValue<Derived> MatrixBase<Derived>::cosh() const
561 {
562  eigen_assert(rows() == cols());
563  typedef typename internal::stem_function<Scalar>::ComplexScalar ComplexScalar;
564  return MatrixFunctionReturnValue<Derived>(derived(), internal::stem_function_cosh<ComplexScalar>);
565 }
566 
567 } // end namespace Eigen
568 
569 #endif // EIGEN_MATRIX_FUNCTION_H
gtsam::internal::rotation
Rot3 rotation(const Pose3 &pose, OptionalJacobian< 3, 6 > H)
Definition: slam/expressions.h:80
Eigen::internal::MatrixFunctionAtomic::compute
MatrixType compute(const MatrixType &A)
Compute matrix function of atomic matrix.
Definition: MatrixFunction.h:64
Eigen::internal::MatrixFunctionAtomic::Scalar
MatrixType::Scalar Scalar
Definition: MatrixFunction.h:34
Eigen::internal::matrix_function_compute_permutation
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:203
Eigen
Namespace containing all symbols from the Eigen library.
Definition: jet.h:637
Eigen::internal::matrix_function_compute< MatrixType, 0 >::run
static void run(const MatA &A, AtomicType &atomic, ResultType &result)
Definition: MatrixFunction.h:393
s
RealScalar s
Definition: level1_cplx_impl.h:126
e
Array< double, 1, 3 > e(1./3., 0.5, 2.)
MatrixType
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
ceres::sin
Jet< T, N > sin(const Jet< T, N > &f)
Definition: jet.h:439
Eigen::MatrixFunctionReturnValue::MatrixFunctionReturnValue
MatrixFunctionReturnValue(const Derived &A, StemFunction f)
Constructor.
Definition: MatrixFunction.h:487
mu
double mu
Definition: testBoundingConstraint.cpp:37
Eigen::MatrixFunctionReturnValue::cols
Index cols() const
Definition: MatrixFunction.h:509
eigen_assert
#define eigen_assert(x)
Definition: Macros.h:1037
Eigen::MatrixFunctionReturnValue::DerivedNested
internal::ref_selector< Derived >::type DerivedNested
Definition: MatrixFunction.h:478
Eigen::internal::matrix_function_compute_block_start
void matrix_function_compute_block_start(const VectorType &clusterSize, VectorType &blockStart)
Compute start of each block using clusterSize.
Definition: MatrixFunction.h:176
Eigen::internal::MatrixFunctionAtomic::m_f
StemFunction * m_f
Definition: MatrixFunction.h:49
B
Definition: test_numpy_dtypes.cpp:299
Eigen::Success
@ Success
Definition: Constants.h:442
Eigen::internal::stem_function
Definition: ForwardDeclarations.h:313
type
Definition: pytypes.h:1525
T
Eigen::Triplet< double > T
Definition: Tutorial_sparse_example.cpp:6
Eigen::internal::matrix_function_find_cluster
ListOfClusters::iterator matrix_function_find_cluster(Index key, ListOfClusters &clusters)
Find cluster in clusters containing some value.
Definition: MatrixFunction.h:108
schurOfA
cout<< "Here is a random 4x4 matrix, A:"<< endl<< A<< endl<< endl;ComplexSchur< MatrixXcf > schurOfA(A, false)
X
#define X
Definition: icosphere.cpp:20
Eigen::MatrixBase::matrixFunction
const MatrixFunctionReturnValue< Derived > matrixFunction(StemFunction f) const
Helper function for the unsupported MatrixFunctions module.
Definition: MatrixFunction.h:529
Eigen::JacobiRotation
Rotation given by a cosine-sine pair.
Definition: ForwardDeclarations.h:283
result
Values result
Definition: OdometryOptimize.cpp:8
Eigen::internal::matrix_function_separation
static const float matrix_function_separation
Maximum distance allowed between eigenvalues to be considered "close".
Definition: MatrixFunction.h:21
rows
int rows
Definition: Tutorial_commainit_02.cpp:1
IsComplex
@ IsComplex
Definition: gtsam/3rdparty/Eigen/blas/common.h:98
ceres::cos
Jet< T, N > cos(const Jet< T, N > &f)
Definition: jet.h:426
Eigen::MatrixFunctionReturnValue::m_f
StemFunction * m_f
Definition: MatrixFunction.h:513
cosh
const EIGEN_DEVICE_FUNC CoshReturnType cosh() const
Definition: ArrayCwiseUnaryOps.h:353
n
int n
Definition: BiCGSTAB_simple.cpp:1
Eigen::internal::traits< MatrixFunctionReturnValue< Derived > >::ReturnType
Derived::PlainObject ReturnType
Definition: MatrixFunction.h:520
Eigen::internal::true_type
Definition: Meta.h:96
A
Definition: test_numpy_dtypes.cpp:298
Eigen::internal::matrix_function_compute_block_atomic
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:242
gtsam.examples.PlanarManipulatorExample.delta
def delta(g0, g1)
Definition: PlanarManipulatorExample.py:45
j
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2
Eigen::internal::matrix_function_compute_map
void matrix_function_compute_map(const EivalsType &eivals, const ListOfClusters &clusters, VectorType &eivalToCluster)
Compute mapping of eigenvalue indices to cluster indices.
Definition: MatrixFunction.h:187
Eigen::MatrixFunctionReturnValue::rows
Index rows() const
Definition: MatrixFunction.h:508
eivals
VectorXcd eivals
Definition: MatrixBase_eigenvalues.cpp:2
Eigen::internal::MatrixFunctionAtomic::MatrixFunctionAtomic
MatrixFunctionAtomic(StemFunction f)
Constructor.
Definition: MatrixFunction.h:40
l
static const Line3 l(Rot3(), 1, 1)
Eigen::PlainObjectBase::rows
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:143
std::swap
void swap(GeographicLib::NearestNeighbor< dist_t, pos_t, distfun_t > &a, GeographicLib::NearestNeighbor< dist_t, pos_t, distfun_t > &b)
Definition: NearestNeighbor.hpp:827
m
Matrix3f m
Definition: AngleAxis_mimic_euler.cpp:1
Eigen::Triplet
A small structure to hold a non zero as a triplet (i,j,value).
Definition: SparseUtil.h:162
Eigen::MatrixFunctionReturnValue::evalTo
void evalTo(ResultType &result) const
Compute the matrix function.
Definition: MatrixFunction.h:494
gtsam::symbol_shorthand::F
Key F(std::uint64_t j)
Definition: inference/Symbol.h:153
key
const gtsam::Symbol key('X', 0)
tree::f
Point2(* f)(const Point3 &, OptionalJacobian< 2, 3 >)
Definition: testExpression.cpp:218
Eigen::internal::MatrixFunctionAtomic
Helper class for computing matrix functions of atomic matrices.
Definition: MatrixFunction.h:30
RealScalar
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:47
Eigen::internal::stem_function::ComplexScalar
std::complex< typename NumTraits< Scalar >::Real > ComplexScalar
Definition: ForwardDeclarations.h:315
Eigen::internal::matrix_function_partition_eigenvalues
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:131
Eigen::internal::matrix_function_compute::run
static void run(const MatrixType &A, AtomicType &atomic, ResultType &result)
Compute the matrix function.
Eigen::internal::MatrixFunctionAtomic::StemFunction
stem_function< Scalar >::type StemFunction
Definition: MatrixFunction.h:35
Eigen::internal::matrix_function_compute_mu
NumTraits< typename MatrixType::Scalar >::Real matrix_function_compute_mu(const MatrixType &A)
Definition: MatrixFunction.h:53
C
Matrix< Scalar, Dynamic, Dynamic > C
Definition: bench_gemm.cpp:50
Eigen::internal::matrix_function_solve_triangular_sylvester
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:274
Eigen::MatrixFunctionReturnValue
Proxy for the matrix function of some matrix (expression).
Definition: ForwardDeclarations.h:305
Eigen::internal::traits< MatrixType >
Eigen::PlainObjectBase::cols
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:145
StemFunction.h
Real
mp::number< mp::cpp_dec_float< 100 >, mp::et_on > Real
Definition: boostmultiprec.cpp:78
Eigen::MatrixFunctionReturnValue::Scalar
Derived::Scalar Scalar
Definition: MatrixFunction.h:474
Eigen::internal::matrix_function_compute_above_diagonal
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:323
Eigen::MatrixFunctionReturnValue::m_A
const DerivedNested m_A
Definition: MatrixFunction.h:512
P
static double P[]
Definition: ellpe.c:68
Eigen::internal::matrix_function_compute
Class for computing matrix functions.
Definition: MatrixFunction.h:367
U
@ U
Definition: testDecisionTree.cpp:342
Eigen::Matrix
The matrix class, also used for vectors and row-vectors.
Definition: 3rdparty/Eigen/Eigen/src/Core/Matrix.h:178
abs
#define abs(x)
Definition: datatypes.h:17
N
#define N
Definition: igam.h:9
internal
Definition: BandTriangularSolver.h:13
VectorType
Definition: FFTW.cpp:65
Eigen::PlainObjectBase< Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > >::setZero
EIGEN_DEVICE_FUNC Matrix< _Scalar, _Rows, _Cols, _Options, _MaxRows, _MaxCols > & setZero(Index size)
Definition: CwiseNullaryOp.h:562
Eigen::internal::matrix_function_compute_cluster_size
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:163
Eigen::internal::matrix_function_permute_schur
void matrix_function_permute_schur(VectorType &permutation, MatrixType &U, MatrixType &T)
Permute Schur decomposition in U and T according to permutation.
Definition: MatrixFunction.h:216
cols
int cols
Definition: Tutorial_commainit_02.cpp:1
Eigen::internal::matrix_function_compute< MatrixType, 1 >::run
static void run(const MatA &A, AtomicType &atomic, ResultType &result)
Definition: MatrixFunction.h:417
Eigen::internal::plain_col_type
Definition: XprHelper.h:614
sinh
const EIGEN_DEVICE_FUNC SinhReturnType sinh() const
Definition: ArrayCwiseUnaryOps.h:339
max
#define max(a, b)
Definition: datatypes.h:20
Eigen::MatrixFunctionReturnValue::StemFunction
internal::stem_function< Scalar >::type StemFunction
Definition: MatrixFunction.h:475
Eigen::ComplexSchur< MatrixType >
Eigen::NumTraits
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:232
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Scalar
SCALAR Scalar
Definition: bench_gemm.cpp:46
Eigen::Index
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
Cols
static const int Cols
Definition: testCallRecord.cpp:32


gtsam
Author(s):
autogenerated on Sun Dec 22 2024 04:12:18