DenseCholesky.h
Go to the documentation of this file.
1 // Copyright (C) 2016-2025 Yixuan Qiu <yixuan.qiu@cos.name>
2 //
3 // This Source Code Form is subject to the terms of the Mozilla
4 // Public License v. 2.0. If a copy of the MPL was not distributed
5 // with this file, You can obtain one at https://mozilla.org/MPL/2.0/.
6 
7 #ifndef SPECTRA_DENSE_CHOLESKY_H
8 #define SPECTRA_DENSE_CHOLESKY_H
9 
10 #include <Eigen/Core>
11 #include <Eigen/Cholesky>
12 #include <stdexcept>
13 
14 #include "../Util/CompInfo.h"
15 
16 namespace Spectra {
17 
33 template <typename Scalar_, int Uplo = Eigen::Lower, int Flags = Eigen::ColMajor>
35 {
36 public:
40  using Scalar = Scalar_;
41 
42 private:
48 
49  const Index m_n;
51  CompInfo m_info; // status of the decomposition
52 
53 public:
62  template <typename Derived>
65  {
66  static_assert(
67  static_cast<int>(Derived::PlainObject::IsRowMajor) == static_cast<int>(Matrix::IsRowMajor),
68  "DenseCholesky: the \"Flags\" template parameter does not match the input matrix (Eigen::ColMajor/Eigen::RowMajor)");
69 
70  if (m_n != mat.cols())
71  throw std::invalid_argument("DenseCholesky: matrix must be square");
72 
73  m_decomp.compute(mat);
74  m_info = (m_decomp.info() == Eigen::Success) ?
77  }
78 
82  Index rows() const { return m_n; }
86  Index cols() const { return m_n; }
87 
92  CompInfo info() const { return m_info; }
93 
100  // y_out = inv(L) * x_in
101  void lower_triangular_solve(const Scalar* x_in, Scalar* y_out) const
102  {
103  MapConstVec x(x_in, m_n);
104  MapVec y(y_out, m_n);
105  y.noalias() = m_decomp.matrixL().solve(x);
106  }
107 
114  // y_out = inv(L') * x_in
115  void upper_triangular_solve(const Scalar* x_in, Scalar* y_out) const
116  {
117  MapConstVec x(x_in, m_n);
118  MapVec y(y_out, m_n);
119  y.noalias() = m_decomp.matrixU().solve(x);
120  }
121 };
122 
123 } // namespace Spectra
124 
125 #endif // SPECTRA_DENSE_CHOLESKY_H
Spectra::DenseCholesky::info
CompInfo info() const
Definition: DenseCholesky.h:92
Spectra::DenseCholesky::m_info
CompInfo m_info
Definition: DenseCholesky.h:51
Spectra::DenseCholesky::upper_triangular_solve
void upper_triangular_solve(const Scalar *x_in, Scalar *y_out) const
Definition: DenseCholesky.h:115
Spectra::CompInfo::Successful
@ Successful
Computation was successful.
Spectra::DenseCholesky< double >::Index
Eigen::Index Index
Definition: DenseCholesky.h:43
x
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition: gnuplot_common_settings.hh:12
Eigen::Success
@ Success
Definition: Constants.h:442
mat
MatrixXf mat
Definition: Tutorial_AdvancedInitialization_CommaTemporary.cpp:1
Spectra::DenseCholesky::m_decomp
Eigen::LLT< Matrix, Uplo > m_decomp
Definition: DenseCholesky.h:50
Spectra::DenseCholesky::rows
Index rows() const
Definition: DenseCholesky.h:82
Spectra::CompInfo
CompInfo
Definition: CompInfo.h:17
Spectra::DenseCholesky
Definition: DenseCholesky.h:34
Eigen::Map
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:94
Spectra::DenseCholesky< double >::Scalar
double Scalar
Definition: DenseCholesky.h:40
y
Scalar * y
Definition: level1_cplx_impl.h:124
Spectra::CompInfo::NumericalIssue
@ NumericalIssue
Eigen::LLT
Standard Cholesky decomposition (LL^T) of a matrix and associated features.
Definition: LLT.h:66
Spectra
Definition: LOBPCGSolver.h:19
Spectra::DenseCholesky::DenseCholesky
DenseCholesky(const Eigen::MatrixBase< Derived > &mat)
Definition: DenseCholesky.h:63
Spectra::DenseCholesky::lower_triangular_solve
void lower_triangular_solve(const Scalar *x_in, Scalar *y_out) const
Definition: DenseCholesky.h:101
Eigen::Matrix
The matrix class, also used for vectors and row-vectors.
Definition: 3rdparty/Eigen/Eigen/src/Core/Matrix.h:178
Eigen::MatrixBase
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Spectra::DenseCholesky::m_n
const Index m_n
Definition: DenseCholesky.h:49
Spectra::DenseCholesky::cols
Index cols() const
Definition: DenseCholesky.h:86
Spectra::CompInfo::NotComputed
@ NotComputed
Eigen::Index
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74


gtsam
Author(s):
autogenerated on Fri Apr 4 2025 03:01:22