Prusa Slicer 2.6.0
Loading...
Searching...
No Matches
Eigen::SparseSolverBase< Derived > Class Template Reference

A base class for sparse solvers. More...

#include <src/eigen/Eigen/src/SparseCore/SparseSolverBase.h>

Inherits Eigen::internal::noncopyable.

Inherited by Eigen::CholmodBase< _MatrixType, Lower, CholmodDecomposition< _MatrixType, Lower > >, Eigen::CholmodBase< _MatrixType, Lower, CholmodSimplicialLDLT< _MatrixType, Lower > >, Eigen::CholmodBase< _MatrixType, Lower, CholmodSimplicialLLT< _MatrixType, Lower > >, Eigen::CholmodBase< _MatrixType, Lower, CholmodSupernodalLLT< _MatrixType, Lower > >, Eigen::IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >, Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >, Eigen::IterativeSolverBase< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >, Eigen::PardisoImpl< PardisoLDLT< MatrixType, Options > >, Eigen::PardisoImpl< PardisoLLT< MatrixType, _UpLo > >, Eigen::PardisoImpl< PardisoLU< MatrixType > >, Eigen::PastixBase< PastixLDLT< _MatrixType, _UpLo > >, Eigen::PastixBase< PastixLLT< _MatrixType, _UpLo > >, Eigen::PastixBase< PastixLU< _MatrixType > >, Eigen::SimplicialCholeskyBase< SimplicialCholesky< _MatrixType, _UpLo, _Ordering > >, Eigen::SimplicialCholeskyBase< SimplicialLDLT< _MatrixType, _UpLo, _Ordering > >, Eigen::SimplicialCholeskyBase< SimplicialLLT< _MatrixType, _UpLo, _Ordering > >, Eigen::SimplicialCholeskyBase< SimplicialLLT< Eigen::SparseMatrix< double >, _UpLo, _Ordering > >, Eigen::SimplicialCholeskyBase< SimplicialLDLT< Eigen::SparseMatrix< double >, _UpLo, _Ordering > >, Eigen::SimplicialCholeskyBase< SimplicialCholesky< Eigen::SparseMatrix< double >, _UpLo, _Ordering > >, Eigen::SimplicialCholeskyBase< SimplicialLDLT< Eigen::SparseMatrix< typename DerivedV::Scalar >, _UpLo, _Ordering > >, Eigen::SimplicialCholeskyBase< SimplicialLLT< Eigen::SparseMatrix< T >, _UpLo, _Ordering > >, Eigen::SimplicialCholeskyBase< SimplicialLDLT< Eigen::SparseMatrix< T >, _UpLo, _Ordering > >, Eigen::SparseLU< Eigen::SparseMatrix< double, Eigen::ColMajor >, Eigen::COLAMDOrdering< int > >, Eigen::SparseLU< Eigen::SparseMatrix< T, Eigen::ColMajor >, Eigen::COLAMDOrdering< int > >, Eigen::SparseQR< Eigen::SparseMatrix< double >, Eigen::COLAMDOrdering< int > >, Eigen::SparseQR< Eigen::SparseMatrix< T >, Eigen::COLAMDOrdering< int > >, Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >, Eigen::CholmodBase< _MatrixType, _UpLo, Derived >, Eigen::IncompleteCholesky< Scalar, _UpLo, _OrderingType >, Eigen::IncompleteLUT< _Scalar, _StorageIndex >, Eigen::IterativeSolverBase< Derived >, Eigen::PardisoImpl< Derived >, Eigen::PastixBase< Derived >, Eigen::SimplicialCholeskyBase< Derived >, and Eigen::SuperLUBase< _MatrixType, Derived >.

+ Collaboration diagram for Eigen::SparseSolverBase< Derived >:

Public Member Functions

 SparseSolverBase ()
 
 ~SparseSolverBase ()
 
Derived & derived ()
 
const Derived & derived () const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
template<typename Rhs , typename Dest >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Protected Attributes

bool m_isInitialized
 

Detailed Description

template<typename Derived>
class Eigen::SparseSolverBase< Derived >

A base class for sparse solvers.

Template Parameters
Derivedthe actual type of the solver.

Constructor & Destructor Documentation

◆ SparseSolverBase()

template<typename Derived >
Eigen::SparseSolverBase< Derived >::SparseSolverBase ( )
inline

Default constructor

73 : m_isInitialized(false)
74 {}
bool m_isInitialized
Definition SparseSolverBase.h:119

◆ ~SparseSolverBase()

template<typename Derived >
Eigen::SparseSolverBase< Derived >::~SparseSolverBase ( )
inline
77 {}

Member Function Documentation

◆ _solve_impl()

template<typename Derived >
template<typename Rhs , typename Dest >
void Eigen::SparseSolverBase< Derived >::_solve_impl ( const SparseMatrixBase< Rhs > &  b,
SparseMatrixBase< Dest > &  dest 
) const
inline
112 {
113 internal::solve_sparse_through_dense_panels(derived(), b.derived(), dest.derived());
114 }
Derived & derived()
Definition SparseSolverBase.h:79
enable_if< Rhs::ColsAtCompileTime!=1 &&Dest::ColsAtCompileTime!=1 >::type solve_sparse_through_dense_panels(const Decomposition &dec, const Rhs &rhs, Dest &dest)
Definition SparseSolverBase.h:23

References Eigen::SparseSolverBase< Derived >::derived(), Eigen::SparseMatrixBase< Derived >::derived(), and Eigen::internal::solve_sparse_through_dense_panels().

+ Here is the call graph for this function:

◆ derived() [1/2]

template<typename Derived >
Derived & Eigen::SparseSolverBase< Derived >::derived ( )
inline
79{ return *static_cast<Derived*>(this); }

Referenced by Eigen::SparseSolverBase< Derived >::_solve_impl(), and Eigen::SparseSolverBase< Derived >::solve().

+ Here is the caller graph for this function:

◆ derived() [2/2]

template<typename Derived >
const Derived & Eigen::SparseSolverBase< Derived >::derived ( ) const
inline
80{ return *static_cast<const Derived*>(this); }

◆ solve() [1/2]

template<typename Derived >
template<typename Rhs >
const Solve< Derived, Rhs > Eigen::SparseSolverBase< Derived >::solve ( const MatrixBase< Rhs > &  b) const
inline
Returns
an expression of the solution x of $ A x = b $ using the current decomposition of A.
See also
compute()
89 {
90 eigen_assert(m_isInitialized && "Solver is not initialized.");
91 eigen_assert(derived().rows()==b.rows() && "solve(): invalid number of rows of the right hand side matrix b");
92 return Solve<Derived, Rhs>(derived(), b.derived());
93 }
#define eigen_assert(x)
Definition Macros.h:579
size_t rows(const T &raster)
Definition MarchingSquares.hpp:55

References Eigen::SparseSolverBase< Derived >::derived(), eigen_assert, and Eigen::SparseSolverBase< Derived >::m_isInitialized.

Referenced by igl::embree::bone_heat(), igl::Frame_field_deformer::compute_optimal_positions(), igl::eigs(), igl::copyleft::comiso::FrameInterpolator::interpolateSymmetric(), igl::slim::solve_weighted_arap(), and igl::copyleft::comiso::NRosyField::solveNoRoundings().

+ Here is the call graph for this function:
+ Here is the caller graph for this function:

◆ solve() [2/2]

template<typename Derived >
template<typename Rhs >
const Solve< Derived, Rhs > Eigen::SparseSolverBase< Derived >::solve ( const SparseMatrixBase< Rhs > &  b) const
inline
Returns
an expression of the solution x of $ A x = b $ using the current decomposition of A.
See also
compute()
102 {
103 eigen_assert(m_isInitialized && "Solver is not initialized.");
104 eigen_assert(derived().rows()==b.rows() && "solve(): invalid number of rows of the right hand side matrix b");
105 return Solve<Derived, Rhs>(derived(), b.derived());
106 }

Member Data Documentation

◆ m_isInitialized

template<typename Derived >
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprotected

The documentation for this class was generated from the following file: