an efficient C++ finite element environment
|
|
Go to the documentation of this file. 1 #ifndef _RHEOLEF_SOLVER_CHOLMOD_H
2 #define _RHEOLEF_SOLVER_CHOLMOD_H
25 #include "rheolef/solver.h"
27 #if defined(_RHEOLEF_HAVE_CHOLMOD)
28 #include <suitesparse/cholmod.h>
29 #pragma GCC diagnostic push
30 #pragma GCC diagnostic ignored "-Weffc++"
31 #include <Eigen/Sparse>
32 #include <Eigen/src/CholmodSupport/CholmodSupport.h>
33 #pragma GCC diagnostic pop
39 template<
class T,
class M>
65 Eigen::CholmodDecomposition<Eigen::SparseMatrix<double> >
72 template<
class T,
class M>
82 template<
class T,
class M>
93 template <
class T,
class M>
99 return new_macro (rep(*
this));
103 #endif // _RHEOLEF_HAVE_CHOLMOD
104 #endif // _RHEOLEF_SOLVER_CHOLMOD_H
Eigen::CholmodDecomposition< Eigen::SparseMatrix< double > > _llt_a
base::determinant_type determinant_type
see the vec page for the full documentation
solver_abstract_rep< T, M > base
vec< T, M > solve(const vec< T, M > &rhs) const
see the csr page for the full documentation
vec< T, M > trans_solve(const vec< T, M > &rhs) const
solver_abstract_rep< T, M > * clone() const
This file is part of Rheolef.
csr< T, M >::size_type size_type
base::size_type size_type
void update_values(const csr< T, M > &a)
determinant_type det() const
solver_cholmod_rep(const csr< T, M > &a, const solver_option &opt=solver_option())
see the solver_option page for the full documentation