Dune Core Modules (2.6.0)

Dune::SPQR< BCRSMatrix< FieldMatrix< T, n, m >, A > > Class Template Referenceabstract

The SPQR direct sparse solver for matrices of type BCRSMatrix. More...

#include <dune/istl/spqr.hh>

Public Types

typedef Dune::BCRSMatrix< FieldMatrix< T, n, m >, A > Matrix
 The matrix type.
 
typedef Dune::ColCompMatrix< MatrixSPQRMatrix
 The corresponding SuperLU Matrix type.
 
typedef ColCompMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > MatrixInitializer
 Type of an associated initializer class.
 
typedef Dune::BlockVector< FieldVector< T, m >, typename A::template rebind< FieldVector< T, m > >::other > domain_type
 The type of the domain of the solver.
 
typedef Dune::BlockVector< FieldVector< T, n >, typename A::template rebind< FieldVector< T, n > >::other > range_type
 The type of the range of the solver.
 
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 SimdScalar< real_typescalar_real_type
 scalar type underlying the field_type
 

Public Member Functions

virtual SolverCategory::Category category () const
 Category of the solver (see SolverCategory::Category)
 
 SPQR (const Matrix &matrix, int verbose=0)
 Construct a solver object from a BCRSMatrix. More...
 
 SPQR (const Matrix &matrix, int verbose, bool)
 Constructor for compatibility with SuperLU standard constructor. More...
 
 SPQR ()
 Default constructor.
 
virtual ~SPQR ()
 Destructor.
 
virtual void apply (domain_type &x, range_type &b, InverseOperatorResult &res)
 Apply inverse operator,. More...
 
virtual void apply (domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
 apply inverse operator, with given convergence criteria. More...
 
void setMatrix (const Matrix &matrix)
 Initialize data from given matrix.
 
void setVerbosity (int v)
 Sets the verbosity level for the solver. More...
 
SuiteSparseQR_factorization< T > * getFactorization ()
 Return the matrix factorization. More...
 
SPQRMatrixgetInternalMatrix ()
 Return the column coppressed matrix. More...
 
void free ()
 Free allocated space. More...
 
const char * name ()
 Get method name.
 
virtual void apply (BlockVector< FieldVector< T, m >, A::template rebind< FieldVector< T, m > >::other > &x, BlockVector< FieldVector< T, n >, A::template rebind< FieldVector< T, n > >::other > &b, InverseOperatorResult &res)=0
 Apply inverse operator,. More...
 
virtual void apply (BlockVector< FieldVector< T, m >, A::template rebind< FieldVector< T, m > >::other > &x, BlockVector< FieldVector< T, n >, A::template rebind< FieldVector< T, n > >::other > &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<typename T, typename A, int n, int m>
class Dune::SPQR< BCRSMatrix< FieldMatrix< T, n, m >, A > >

The SPQR direct sparse solver for matrices of type BCRSMatrix.

Specialization for the Dune::BCRSMatrix. SPQR will always go double precision and supports complex numbers too (use std::complex<double> for that).

Template Parameters
TNumber type. Only double and std::complex<double> is supported
ASTL-compatible allocator type
nNumber of rows in a matrix block
mNumber of columns in a matrix block
Note
This will only work if dune-istl has been configured to use SPQR

Member Function Documentation

◆ apply() [1/2]

virtual void Dune::InverseOperator< BlockVector< FieldVector< T, m >, A::template rebind< FieldVector< T, m > >::other > , BlockVector< FieldVector< T, n >, A::template rebind< FieldVector< T, n > >::other > >::apply ( X &  x,
Y &  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]

virtual void Dune::InverseOperator< BlockVector< FieldVector< T, m >, A::template rebind< FieldVector< T, m > >::other > , BlockVector< FieldVector< T, n >, A::template rebind< FieldVector< T, n > >::other > >::apply ( X &  x,
Y &  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

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