KASKADE 7 development version
|
Classes and functions for elimination methods. More...
Classes | |
class | Kaskade::DirectSolver< Domain_, Range_ > |
class | Kaskade::DirectPreconditioner< Op > |
class | Kaskade::Factorization< Scalar > |
Abstract base class for matrix factorizations. More... | |
struct | Kaskade::Creator< Factorization< Scalar > > |
Abstract base class for factorization creators to be plugged into a factorization factory. More... | |
class | Kaskade::MUMPSFactorization< Scalar > |
Factorization of sparse linear systems with mumps. More... | |
Enumerations | |
enum class | DirectType { DirectType::ANY =-1 , DirectType::UMFPACK , DirectType::PARDISO , DirectType::MUMPS , DirectType::SUPERLU , DirectType::UMFPACK3264 , DirectType::UMFPACK64 } |
Available direct solvers for linear equation systems. More... | |
Functions | |
template<class Scalar , class Index > | |
std::unique_ptr< Factorization< Scalar > > | Kaskade::getFactorization (DirectType directType, MatrixAsTriplet< Scalar, Index > const &A, FactorizationOptions options) |
Creates a factorization of the given triplet matrix. More... | |
template<class Scalar > | |
std::unique_ptr< Factorization< Scalar > > | Kaskade::getFactorization (DirectType directType, Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, 1, 1 > > const &A, FactorizationOptions options) |
Creates a factorization of the given BCRSmatrix. More... | |
template<class Scalar , int n, int m, class Index > | |
std::unique_ptr< Factorization< Scalar > > | Kaskade::getFactorization (DirectType directType, NumaBCRSMatrix< Dune::FieldMatrix< Scalar, n, m >, Index > const &A, FactorizationOptions options) |
Creates a factorization of the given NumaBCRS matrix. More... | |
Classes and functions for elimination methods.
Kaskade 7 provides a dynamic abstract interface to direct solvers for linear equation systems \( Ax = b \). The interface and available direct solvers are described here.
Different solvers have been integrated by providing glue code registering them dynamically with the solver factory. Only solvers that are actually installed are available. This can be queried at run time at the factory.
Direct type: UMFPACK, UMFPACK64
UMFPack is a general-purpose direct solver computing an LU decomposition with pivoting, included in the SuiteSparse software. It is not the fastest solver around, but robust and easy to work with. Unless performance or memory consumption are critical, it is a good choice to start with.
In 3D diffusion problems, UMFPack works well up to more than 100k degrees of freedom. Standard UMFPACK works with 32 bit int indices, and can therefore fail on very large problems. In that case, use UMFPACK64, which works with 64 bit long indices.
UMFPack supports the following options that can be set in the FactorizationOptions object passed when constructing a solver:
refinementIterations
: number of refinement iterations to perform. Default is 2. With at least one refinement iteration, LU decomposition with pivoting is numerically stable. In many cases, e.g. in preconditioners, with well-behaved symmetric positive definite systems, or if discretization error is anyways larger, setting this to 0 can reduce the computational cost of the solution phase (not the factorization, of course) by a factor of 3 compared to the default, without significant loss of accuracy.
|
strong |
Available direct solvers for linear equation systems.
Not all of these need be actually available, that depends on which solvers are included in the installation. But these are the lookup keys for creating solvers via the factorization factory.
std::unique_ptr< Factorization< Scalar > > Kaskade::getFactorization | ( | DirectType | directType, |
Dune::BCRSMatrix< Dune::FieldMatrix< Scalar, 1, 1 > > const & | A, | ||
FactorizationOptions | options | ||
) |
Creates a factorization of the given BCRSmatrix.
Definition at line 217 of file factorization.hh.
std::unique_ptr< Factorization< Scalar > > Kaskade::getFactorization | ( | DirectType | directType, |
MatrixAsTriplet< Scalar, Index > const & | A, | ||
FactorizationOptions | options | ||
) |
Creates a factorization of the given triplet matrix.
Definition at line 192 of file factorization.hh.
Referenced by Kaskade::getFactorization(), and Kaskade::PartialDirectPreconditioner< Op >::PartialDirectPreconditioner().
std::unique_ptr< Factorization< Scalar > > Kaskade::getFactorization | ( | DirectType | directType, |
NumaBCRSMatrix< Dune::FieldMatrix< Scalar, n, m >, Index > const & | A, | ||
FactorizationOptions | options | ||
) |
Creates a factorization of the given NumaBCRS matrix.
Definition at line 231 of file factorization.hh.