Dune Core Modules (2.4.1)
Iterative Solvers supporting block recursive matrix and vector classes at compile time. More...
Modules | |
Sparse Matrix and Vector classes | |
Matrix and Vector classes that support a block recursive structure capable of representing the natural structure from Finite Element discretisations. | |
Iterative Solvers | |
Files | |
file | superlu.hh |
Classes for using SuperLU with ISTL matrices. | |
file | umfpack.hh |
Classes for using UMFPack with ISTL matrices. | |
Classes | |
class | Dune::ILUSubdomainSolver< M, X, Y > |
base class encapsulating common algorithms of ILU0SubdomainSolver and ILUNSubdomainSolver. More... | |
class | Dune::ILU0SubdomainSolver< M, X, Y > |
Exact subdomain solver using ILU(p) with appropriate p. More... | |
class | Dune::ISTLError |
derive error class from the base class in common More... | |
class | Dune::BCRSMatrixError |
Error specific to BCRSMatrix. More... | |
class | Dune::ImplicitModeOverflowExhausted |
The overflow error used during implicit BCRSMatrix construction was exhausted. More... | |
class | Dune::SeqOverlappingSchwarz< M, X, TM, TD, TA > |
Sequential overlapping Schwarz preconditioner. More... | |
class | Dune::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > > |
SuperLu Solver. More... | |
class | Dune::UMFPack< Matrix > |
Use the UMFPack package to directly solve linear systems – empty default class. More... | |
class | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > > |
The UMFPack direct sparse solver for matrices of type BCRSMatrix. More... | |
Typedefs | |
typedef Dune::BCRSMatrix< FieldMatrix< T, n, m >, A > | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::Matrix |
The matrix type. | |
typedef Dune::ColCompMatrix< Matrix > | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::UMFPackMatrix |
The corresponding SuperLU Matrix type. | |
typedef ColCompMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > | Dune::UMFPack< 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 > | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::domain_type |
The type of the domain of the solver. | |
typedef Dune::BlockVector< FieldVector< T, n >, typename A::template rebind< FieldVector< T, n > >::other > | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::range_type |
The type of the range of the solver. | |
Functions | |
Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::UMFPack (const Matrix &matrix, int verbose=0) | |
Construct a solver object from a BCRSMatrix. More... | |
Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::UMFPack (const Matrix &matrix, int verbose, bool) | |
Constructor for compatibility with SuperLU standard constructor. More... | |
Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::UMFPack () | |
default constructor | |
Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::UMFPack (const Matrix &mat_, const char *file, int verbose=0) | |
Try loading a decomposition from file and do a decomposition if unsuccessful. More... | |
Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::UMFPack (const char *file, int verbose=0) | |
try loading a decomposition from file More... | |
virtual void | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::apply (domain_type &x, range_type &b, InverseOperatorResult &res) |
Apply inverse operator,. More... | |
virtual void | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::apply (domain_type &x, range_type &b, double reduction, InverseOperatorResult &res) |
apply inverse operator, with given convergence criteria. More... | |
void | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::apply (T *x, T *b) |
additional apply method with c-arrays in analogy to superlu More... | |
void | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::setOption (unsigned int option, double value) |
Set UMFPack-specific options. More... | |
void | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::saveDecomposition (const char *file) |
saves a decomposition to a file More... | |
void | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::setMatrix (const Matrix &matrix) |
Initialize data from given matrix. | |
void | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::setVerbosity (int v) |
sets the verbosity level for the UMFPack solver More... | |
void | Dune::UMFPack< BCRSMatrix< FieldMatrix< T, n, m >, A > >::free () |
free allocated space. More... | |
template<class S > | |
std::size_t | Dune::ILUSubdomainSolver< M, X, Y >::copyToLocalMatrix (const M &A, S &rowset) |
Copy the local part of the global matrix to ILU. More... | |
template<class S > | |
void | Dune::ILU0SubdomainSolver< M, X, Y >::setSubMatrix (const M &A, S &rowset) |
Set the data of the local problem. More... | |
template<class S > | |
void | Dune::ILUNSubdomainSolver< M, X, Y >::setSubMatrix (const M &A, S &rowset) |
Set the data of the local problem. More... | |
void | Dune::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > >::free () |
free allocated space. More... | |
Dune::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > >::SuperLU (const Matrix &mat, bool verbose=false, bool reusevector=true) | |
Constructs the SuperLU solver. More... | |
Dune::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > >::SuperLU () | |
Empty default constructor. More... | |
void | Dune::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > >::setMatrix (const Matrix &mat) |
Initialize data from given matrix. | |
void | Dune::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > >::apply (domain_type &x, range_type &b, InverseOperatorResult &res) |
Apply inverse operator,. More... | |
void | Dune::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > >::apply (T *x, T *b) |
Apply SuperLu to C arrays. | |
Detailed Description
Iterative Solvers supporting block recursive matrix and vector classes at compile time.
The Iterative Solver Template Library applies generic programming in C++ to the domain of iterative solvers of linear systems stemming from finite element discretizations. Those discretizations exhibit a lot of structure, e.g:
- Certain discretizations for systems of PDEs or higher order methods result in matrices where individual entries are replaced by small blocks, say of size \(2\times 2\) or \(4\times 4\). Dense blocks of different sizes e.g. arise in \(hp\) Discontinuous Galerkin discretization methods. It is straightforward and efficient to treat these small dense blocks as fully coupled and solve them with direct methods within the iterative method.
- Equation-wise ordering for systems results in matrices having an \(n\times n\) block structure where \(n\) corresponds to the number of variables in the PDE and the blocks themselves are large and sparse. As an example we mention the Stokes system.
- Other discretisation, e.~g. those of reaction/diffusion systems, produce sparse matrices whose blocks are sparse matrices of small dense blocks,
Our matrix and vector interface supports a block recursive structure. Each sparse matrix entry can itself be either a sparse or a small dense matrix.
The solvers use this recursive block structure via template meta programming at compile time.
ISTL consists of the matrix and vector API and the solvers which use the Preconditioners preconditioners.
Function Documentation
◆ apply() [1/4]
|
inlinevirtual |
apply inverse operator, with given convergence criteria.
- Warning
- Right hand side b may be overwritten!
- Parameters
-
x The left hand side to store the result in. b The right hand side reduction The minimum defect reduction to achieve. res Object to store the statistics about applying the operator.
References DUNE_UNUSED_PARAMETER.
◆ apply() [2/4]
void Dune::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > >::apply | ( | domain_type & | x, |
range_type & | b, | ||
InverseOperatorResult & | res | ||
) |
Apply inverse operator,.
- Warning
- Note: right hand side b may be overwritten!
- Parameters
-
x The left hand side to store the result in. b The right hand side res Object to store the statistics about applying the operator.
◆ apply() [3/4]
|
inlinevirtual |
Apply inverse operator,.
- Warning
- Note: right hand side b may be overwritten!
- Parameters
-
x The left hand side to store the result in. b The right hand side res Object to store the statistics about applying the operator.
References Dune::InverseOperatorResult::converged, Dune::block_vector_unmanaged< B, A >::dim(), DUNE_THROW, Dune::InverseOperatorResult::elapsed, and Dune::InverseOperatorResult::iterations.
◆ apply() [4/4]
|
inline |
additional apply method with c-arrays in analogy to superlu
- Parameters
-
x solution array b rhs array
◆ copyToLocalMatrix()
|
protected |
Copy the local part of the global matrix to ILU.
- Parameters
-
A The global matrix. rowset The global indices of the local problem.
◆ free() [1/2]
void Dune::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > >::free |
free allocated space.
- Warning
- later calling apply will result in an error.
◆ free() [2/2]
|
inline |
free allocated space.
- Warning
- later calling apply will result in an error.
◆ saveDecomposition()
|
inline |
◆ setOption()
|
inline |
Set UMFPack-specific options.
This method allows to set various options that control the UMFPack solver. More specifically, it allows to set values in the UMF_Control array. Please see the UMFPack documentation for a list of possible options and values.
- Parameters
-
option Entry in the UMF_Control array, e.g., UMFPACK_IRSTEP value Corresponding value
- Exceptions
-
RangeError If nonexisting option was requested
References DUNE_THROW.
◆ setSubMatrix() [1/2]
void Dune::ILU0SubdomainSolver< M, X, Y >::setSubMatrix | ( | const M & | A, |
S & | rowset | ||
) |
Set the data of the local problem.
- Parameters
-
A The global matrix. rowset The global indices of the local problem.
- Template Parameters
-
S The type of the set with the indices.
◆ setSubMatrix() [2/2]
void Dune::ILUNSubdomainSolver< M, X, Y >::setSubMatrix | ( | const M & | A, |
S & | rowset | ||
) |
Set the data of the local problem.
- Parameters
-
A The global matrix. rowset The global indices of the local problem.
- Template Parameters
-
S The type of the set with the indices.
◆ setVerbosity()
|
inline |
sets the verbosity level for the UMFPack solver
- Parameters
-
v verbosity level The following levels are implemented: 0 - only error messages 1 - a bit of statistics on decomposition and solution 2 - lots of statistics on decomposition and solution
◆ SuperLU() [1/2]
Dune::SuperLU< BCRSMatrix< FieldMatrix< T, n, m >, A > >::SuperLU |
Empty default constructor.
Use setMatrix to tell SuperLU for what matrix it solves. Using this constructor no vectors will be reused.
◆ SuperLU() [2/2]
|
explicit |
Constructs the SuperLU solver.
During the construction the matrix will be decomposed. That means that in each apply call forward and backward substitutions take place (and no decomposition).
- Parameters
-
mat The matrix of the system to solve. verbose If true some statistics are printed. reusevector Default value is true. If true the two vectors are allocate in the first call to apply. These get resused in subsequent calls to apply and are deallocated in the destructor. If false these vectors are allocated at the beginning and deallocated at the end of each apply method. This allows using the same instance of superlu from different threads.
◆ UMFPack() [1/4]
|
inline |
try loading a decomposition from file
- Parameters
-
file the decomposition file verbose the verbosity level
- Exceptions
-
Dune::Exception When not being able to load the file. Does not need knowledge of the actual matrix!
References DUNE_THROW.
◆ UMFPack() [2/4]
|
inline |
Try loading a decomposition from file and do a decomposition if unsuccessful.
- Parameters
-
mat_ the matrix to decompose when no decoposition file found file the decomposition file verbose the verbosity level
Use saveDecomposition(char* file) for manually storing a decomposition. This constructor will decompose mat_ and store the result to file if no file wasn't found in the first place. Thus, if you always use this you will only compute the decomposition once (and when you manually deleted the decomposition file).
◆ UMFPack() [3/4]
|
inline |
Constructor for compatibility with SuperLU standard constructor.
This computes the matrix decomposition, and may take a long time (and use a lot of memory).
- Parameters
-
mat_ the matrix to solve for verbose [0..2] set the verbosity level, defaults to 0
◆ UMFPack() [4/4]
|
inline |
Construct a solver object from a BCRSMatrix.
This computes the matrix decomposition, and may take a long time (and use a lot of memory).
- Parameters
-
mat_ the matrix to solve for verbose [0..2] set the verbosity level, defaults to 0