Dune Core Modules (unstable)

Dune::Cholmod< Vector, Index > Class Template Referenceabstract

Dune wrapper for SuiteSparse/CHOLMOD solver. More...

#include <dune/istl/cholmod.hh>

Public Types

typedef Vector domain_type
 Type of the domain of the operator to be inverted.
 
typedef Vector range_type
 Type of the range of the operator to be inverted.
 
typedef X::field_type field_type
 The field type of the operator.
 
typedef FieldTraits< field_type >::real_type real_type
 The real type of the field type (is the same if using real numbers, but differs for std::complex)
 
typedef Simd::Scalar< real_typescalar_real_type
 scalar type underlying the field_type
 

Public Member Functions

 Cholmod ()
 Default constructor. More...
 
 ~Cholmod ()
 Destructor. More...
 
void apply (Vector &x, Vector &b, [[maybe_unused]] double reduction, InverseOperatorResult &res)
 simple forward to apply(X&, Y&, InverseOperatorResult&)
 
void apply (Vector &x, Vector &b, InverseOperatorResult &res)
 solve the linear system Ax=b (possibly with respect to some ignore field) More...
 
template<class Matrix >
void setMatrix (const Matrix &matrix)
 Set matrix without ignore nodes. More...
 
template<class Matrix , class Ignore >
void setMatrix (const Matrix &matrix, const Ignore *ignore)
 Set matrix and ignore nodes. More...
 
virtual SolverCategory::Category category () const
 Category of the solver (see SolverCategory::Category)
 
cholmod_common & cholmodCommonObject ()
 return a reference to the CHOLMOD common object for advanced option settings More...
 
cholmod_factor & cholmodFactor ()
 The CHOLMOD data structure that stores the factorization. More...
 
const cholmod_factor & cholmodFactor () const
 The CHOLMOD data structure that stores the factorization. More...
 
virtual void apply (Vector &x, Vector &b, double reduction, InverseOperatorResult &res)=0
 apply inverse operator, with given convergence criteria. More...
 

Protected Member Functions

void printHeader (std::ostream &s) const
 helper function for printing header of solver output
 
void printOutput (std::ostream &s, const CountType &iter, const DataType &norm, const DataType &norm_old) const
 helper function for printing solver output
 
void printOutput (std::ostream &s, const CountType &iter, const DataType &norm) const
 helper function for printing solver output
 

Detailed Description

template<class Vector, class Index = int>
class Dune::Cholmod< Vector, Index >

Dune wrapper for SuiteSparse/CHOLMOD solver.

This class implements an InverseOperator between Vector types

Template Parameters
VectorData type for solution and right-hand-side vectors
IndexType used by CHOLMOD for indices. Must be either 'int' or 'SuiteSparse_long' (which is usually long int).

Constructor & Destructor Documentation

◆ Cholmod()

template<class Vector , class Index = int>
Dune::Cholmod< Vector, Index >::Cholmod ( )
inline

Default constructor.

Calls the the cholmod runtime, see CHOLMOD doc

◆ ~Cholmod()

template<class Vector , class Index = int>
Dune::Cholmod< Vector, Index >::~Cholmod ( )
inline

Destructor.

Free space and calls cholmod_finish, see CHOLMOD doc

Member Function Documentation

◆ apply() [1/2]

virtual void Dune::InverseOperator< Vector , Vector >::apply ( Vector &  x,
Vector &  b,
double  reduction,
InverseOperatorResult res 
)
pure virtualinherited

apply inverse operator, with given convergence criteria.

Warning
Right hand side b may be overwritten!
Parameters
xThe left hand side to store the result in.
bThe right hand side
reductionThe minimum defect reduction to achieve.
resObject to store the statistics about applying the operator.
Exceptions
SolverAbortWhen the solver detects a problem and cannot continue

◆ apply() [2/2]

template<class Vector , class Index = int>
void Dune::Cholmod< Vector, Index >::apply ( Vector &  x,
Vector &  b,
InverseOperatorResult res 
)
inlinevirtual

solve the linear system Ax=b (possibly with respect to some ignore field)

The method assumes that setMatrix() was called before In the case of a given ignore field the corresponding entries of both in x and b will stay untouched in this method.

Implements Dune::InverseOperator< Vector, Vector >.

References Dune::InverseOperatorResult::converged, DUNE_THROW, Dune::flatVectorForEach(), Dune::InverseOperatorResult::iterations, and Dune::Hybrid::max.

◆ cholmodCommonObject()

template<class Vector , class Index = int>
cholmod_common& Dune::Cholmod< Vector, Index >::cholmodCommonObject ( )
inline

return a reference to the CHOLMOD common object for advanced option settings

The CHOLMOD common object stores all parameters and options for the solver to run and can be modified in several ways, see CHOLMOD Userguide for further information.

◆ cholmodFactor() [1/2]

template<class Vector , class Index = int>
cholmod_factor& Dune::Cholmod< Vector, Index >::cholmodFactor ( )
inline

The CHOLMOD data structure that stores the factorization.

Access to this is necessary for the more advanced features of CHOLMOD. You need to know what you are doing!

◆ cholmodFactor() [2/2]

template<class Vector , class Index = int>
const cholmod_factor& Dune::Cholmod< Vector, Index >::cholmodFactor ( ) const
inline

The CHOLMOD data structure that stores the factorization.

Access to this is necessary for the more advanced features of CHOLMOD. You need to know what you are doing!

◆ setMatrix() [1/2]

template<class Vector , class Index = int>
template<class Matrix >
void Dune::Cholmod< Vector, Index >::setMatrix ( const Matrix matrix)
inline

Set matrix without ignore nodes.

This method forwards a nullptr to the setMatrix method with ignore nodes

◆ setMatrix() [2/2]

template<class Vector , class Index = int>
template<class Matrix , class Ignore >
void Dune::Cholmod< Vector, Index >::setMatrix ( const Matrix matrix,
const Ignore *  ignore 
)
inline

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]matrixMatrix to be decomposed. In BCRS compatible form
[in]ignorePointer to a compatible BitVector

References Dune::flatMatrixForEach().


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 27, 22:29, 2024)