Dune wrapper for SuiteSparse/CHOLMOD solver.
More...
#include <dune/istl/cholmod.hh>
|
| Cholmod () |
| Default constructor. More...
|
|
| ~Cholmod () |
| Destructor. More...
|
|
| Cholmod (const Cholmod &)=delete |
|
Cholmod & | operator= (const Cholmod &)=delete |
|
void | apply (X &x, B &b, double reduction, InverseOperatorResult &res) |
| simple forward to apply(X&, Y&, InverseOperatorResult&) More...
|
|
void | apply (X &x, B &b, InverseOperatorResult &res) |
| solve the linear system Ax=b (possibly with respect to some ignore field) More...
|
|
void | setMatrix (const BCRSMatrix< FieldMatrix< T, k, k >> &matrix) |
| Set matrix without ignore nodes. More...
|
|
template<class Ignore > |
void | setMatrix (const BCRSMatrix< FieldMatrix< T, k, k >> &matrix, const Ignore *ignore) |
| Set matrix and ignore nodes. More...
|
|
virtual SolverCategory::Category | category () const |
| Category of the solver (see SolverCategory::Category) More...
|
|
|
void | printHeader (std::ostream &s) const |
| helper function for printing header of solver output More...
|
|
void | printOutput (std::ostream &s, const CountType &iter, const DataType &norm, const DataType &norm_old) const |
| helper function for printing solver output More...
|
|
void | printOutput (std::ostream &s, const CountType &iter, const DataType &norm) const |
| helper function for printing solver output More...
|
|
template<class T, class A, int k>
class Dune::Cholmod< BlockVector< FieldVector< T, k >, A > >
Dune wrapper for SuiteSparse/CHOLMOD solver.
This class implements an InverseOperator between BlockVector types
template<class T , class A , int k>
◆ domain_type
Type of the domain of the operator to be inverted.
◆ field_type
The field type of the operator.
◆ range_type
Type of the range of the operator to be inverted.
◆ real_type
The real type of the field type (is the same if using real numbers, but differs for std::complex)
◆ scalar_real_type
scalar type underlying the field_type
template<class T , class A , int k>
◆ anonymous enum
◆ Cholmod() [1/2]
template<class T , class A , int k>
Default constructor.
Calls the the cholmod runtime, see CHOLMOD doc
◆ ~Cholmod()
template<class T , class A , int k>
Destructor.
Free space and calls cholmod_finish, see CHOLMOD doc
◆ Cholmod() [2/2]
template<class T , class A , int k>
◆ apply() [1/2]
template<class T , class A , int k>
◆ apply() [2/2]
template<class T , class A , int k>
◆ category()
template<class T , class A , int k>
◆ operator=()
template<class T , class A , int k>
◆ printHeader()
helper function for printing header of solver output
◆ printOutput() [1/2]
helper function for printing solver output
◆ printOutput() [2/2]
void Dune::InverseOperator< BlockVector< FieldVector< T, k >, A > , BlockVector< FieldVector< T, k >, A > >::printOutput |
( |
std::ostream & |
s, |
|
|
const CountType & |
iter, |
|
|
const DataType & |
norm, |
|
|
const DataType & |
norm_old |
|
) |
| const |
|
inlineprotectedinherited |
helper function for printing solver output
◆ setMatrix() [1/2]
template<class T , class A , int k>
Set matrix without ignore nodes.
This method forwards a nullptr to the setMatrix method with ignore nodes
◆ setMatrix() [2/2]
template<class T , class A , int k>
template<class Ignore >
Set matrix and ignore nodes.
The input matrix is copied to CHOLMOD compatible form. It is possible to ignore some degrees of freedom, provided an ignore field is given with same block structure like the BlockVector template of the class.
The ignore field causes the method to ignore both rows and cols of the matrix and therefore operates only on the reduced quadratic matrix. In case of, e.g., Dirichlet values the user has to take care of proper adjusting of the rhs before calling apply().
Decomposing the matrix at the end takes a lot of time
- Parameters
-
[in] | matrix | Matrix to be decomposed. In BCRS compatible form |
[in] | ignore | Pointer to a compatible BitVector |
The documentation for this class was generated from the following file: