Dune Core Modules (2.6.0)

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

SuperLu Solver. More...

#include <dune/istl/superlu.hh>

Public Types

typedef Dune::BCRSMatrix< FieldMatrix< T, n, m >, A > Matrix
 The matrix type.
 
typedef Dune::SuperLUMatrix< MatrixSuperLUMatrix
 The corresponding SuperLU Matrix type.
 
typedef SuperMatrixInitializer< 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)
 
 SuperLU (const Matrix &mat, bool verbose=false, bool reusevector=true)
 Constructs the SuperLU solver. More...
 
 SuperLU ()
 Empty default constructor. More...
 
void apply (domain_type &x, range_type &b, InverseOperatorResult &res)
 Apply inverse operator,. More...
 
void apply (domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
 apply inverse operator, with given convergence criteria. More...
 
void apply (T *x, T *b)
 Apply SuperLu to C arrays.
 
void setMatrix (const Matrix &mat)
 Initialize data from given matrix.
 
void free ()
 free allocated space. 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, 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::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > >

SuperLu Solver.

Uses the well known SuperLU package to solve the system.

SuperLU supports single and double precision floating point and complex numbers. Unfortunately these cannot be used at the same time. Therfore users must set SUPERLU_NTYPE (0: float, 1: double, 2: std::complex<float>, 3: std::complex<double>) if the numeric type should be different from double.

Member Function Documentation

◆ apply() [1/3]

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/3]

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

◆ apply() [3/3]

template<typename T , typename A , int n, int m>
void Dune::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > >::apply ( domain_type x,
range_type b,
double  reduction,
InverseOperatorResult res 
)
inline

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

References Dune::Std::apply(), and DUNE_UNUSED_PARAMETER.


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 (Mar 28, 23:30, 2024)