Dune Core Modules (2.7.1)

Dune::Cholmod< BlockVector< FieldVector< T, k >, A > > Class Template Referenceabstract

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

#include <dune/istl/cholmod.hh>

Public Types

typedef BlockVector< FieldVector< T, k >, A > domain_type
 Type of the domain of the operator to be inverted.
 
typedef BlockVector< FieldVector< T, k >, A > 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 (X &x, B &b, double reduction, InverseOperatorResult &res)
 simple forward to apply(X&, Y&, InverseOperatorResult&)
 
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)
 
virtual void apply (BlockVector< FieldVector< T, k >, A > &x, BlockVector< FieldVector< T, k >, A > &b, InverseOperatorResult &res)=0
 Apply inverse operator,. More...
 
virtual void apply (BlockVector< FieldVector< T, k >, A > &x, BlockVector< FieldVector< T, k >, A > &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 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

Constructor & Destructor Documentation

◆ Cholmod()

template<class T , class A , int k>
Dune::Cholmod< BlockVector< FieldVector< T, k >, A > >::Cholmod ( )
inline

Default constructor.

Calls the the cholmod runtime, see CHOLMOD doc

◆ ~Cholmod()

template<class T , class A , int k>
Dune::Cholmod< BlockVector< FieldVector< T, k >, A > >::~Cholmod ( )
inline

Destructor.

Free space and calls cholmod_finish, see CHOLMOD doc

Member Function Documentation

◆ apply() [1/3]

virtual void Dune::InverseOperator< BlockVector< FieldVector< T, k >, A > , BlockVector< FieldVector< T, k >, A > >::apply ( BlockVector< FieldVector< T, k >, A > &  x,
BlockVector< FieldVector< T, k >, A > &  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/3]

virtual void Dune::InverseOperator< BlockVector< FieldVector< T, k >, A > , BlockVector< FieldVector< T, k >, A > >::apply ( BlockVector< FieldVector< T, k >, A > &  x,
BlockVector< FieldVector< T, k >, A > &  b,
InverseOperatorResult res 
)
pure virtualinherited

Apply inverse operator,.

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

◆ apply() [3/3]

template<class T , class A , int k>
void Dune::Cholmod< BlockVector< FieldVector< T, k >, A > >::apply ( X x,
B b,
InverseOperatorResult res 
)
inline

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.

References DUNE_THROW.

◆ setMatrix() [1/2]

template<class T , class A , int k>
void Dune::Cholmod< BlockVector< FieldVector< T, k >, A > >::setMatrix ( const BCRSMatrix< FieldMatrix< T, k, k >> &  matrix)
inline

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 >
void Dune::Cholmod< BlockVector< FieldVector< T, k >, A > >::setMatrix ( const BCRSMatrix< FieldMatrix< T, k, k >> &  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

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 29, 22:29, 2024)