version 3.9-dev

Linear solvers and helpers. More...

Description

Files

file  dunevectors.hh
 Helper to extract native Dune vector types from particular Dumux types.
 
file  helmholtzoperator.hh
 Scalar Helmholtz operator to be used as a solver component.
 
file  istlsolverfactorybackend.hh
 Provides a generic linear solver based on the ISTL that chooses the solver and preconditioner at runtime.
 
file  istlsolverregistry.hh
 The specialized Dumux macro and tag for the ISTL registry to choose the solver and preconditioner at runtime.
 
file  istlsolvers.hh
 Linear solvers from dune-istl.
 
file  linearalgebratraits.hh
 Define traits for linear algebra backends.
 
file  linearsolverparameters.hh
 Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL.
 
file  linearsolvertraits.hh
 Define traits for linear solvers.
 
file  matrixconverter.hh
 A helper class that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix.
 
file  parallelhelpers.hh
 Provides a helper class for nonoverlapping decomposition.
 
file  parallelmatrixadapter.hh
 A parallel version of a linear operator.
 
file  linear/pdesolver.hh
 A high-level solver interface for a linear PDE solver.
 
file  preconditioners.hh
 Dumux preconditioners for iterative solvers.
 
file  scalarproducts.hh
 Classes to compute scalar products.
 
file  scotchbackend.hh
 An interface to the scotch library for matrix reordering.
 
file  seqsolverbackend.hh
 Dumux sequential linear solver backends.
 
file  solver.hh
 Base class for linear solvers.
 
file  solvercategory.hh
 Solver category.
 
file  stokes_solver.hh
 Preconditioned iterative solver for the incompressible Stokes problem.
 
file  symmetrize_constraints.hh
 Helper function to symmetrize row-constraints in constrained linear systems.
 

Namespaces

namespace  Dumux::Detail
 Distance implementation details.
 

Classes

class  Dumux::Detail::IstlIterativeLinearSolver< LinearSolverTraits, LinearAlgebraTraits, InverseOperator, PreconditionerFactory, convertMultiTypeLATypes >
 Standard dune-istl iterative linear solvers. More...
 
class  Dumux::Detail::DirectIstlSolver< LSTraits, LATraits, Solver, convertMultiTypeVectorAndMatrix >
 Direct dune-istl linear solvers. More...
 
class  Dumux::LinearSolverParameters< LinearSolverTraits >
 Generates a parameter tree required for the linear solvers and precondioners of the Dune ISTL. More...
 
class  Dumux::MatrixConverter< MultiTypeBlockMatrix, Scalar >
 A helper class that converts a Dune::MultiTypeBlockMatrix into a plain Dune::BCRSMatrix TODO: allow block sizes for BCRSMatrix other than 1x1 ? More...
 
class  Dumux::VectorConverter< MultiTypeBlockVector, Scalar >
 A helper class that converts a Dune::MultiTypeBlockVector into a plain Dune::BlockVector and transfers back values. More...
 
class  Dumux::ParallelMatrixHelper< Matrix, GridView, RowDofMapper, rowDofCodim >
 Helper class for adding up matrix entries for border entities. More...
 
class  Dumux::LinearPDESolver< Assembler, LinearSolver, Comm >
 An implementation of a linear PDE solver. More...
 
class  Dumux::SeqUzawa< M, X, Y, l >
 A preconditioner based on the Uzawa algorithm for saddle-point problems of the form \( \begin{pmatrix} A & B \\ C & D \end{pmatrix} \begin{pmatrix} u\\ p \end{pmatrix} = \begin{pmatrix} f\\ g \end{pmatrix} \). More...
 
class  Dumux::ScotchBackend< IndexType >
 A reordering backend using the scotch library. More...
 
class  Dumux::IterativePreconditionedSolverImpl
 A general solver backend allowing arbitrary preconditioners and solvers. More...
 
class  Dumux::ExplicitDiagonalSolver
 Solver for simple block-diagonal matrices (e.g. from explicit time stepping schemes) More...
 
class  Dumux::UzawaBiCGSTABBackend< LinearSolverTraits >
 A Uzawa preconditioned BiCGSTAB solver for saddle-point problems. More...
 
class  Dumux::BlockDiagILU0Preconditioner< M, X, Y, blockLevel >
 A simple ilu0 block diagonal preconditioner. More...
 
class  Dumux::BlockDiagILU0BiCGSTABSolver
 A simple ilu0 block diagonal preconditioned BiCGSTABSolver. More...
 
class  Dumux::BlockDiagILU0RestartedGMResSolver
 A simple ilu0 block diagonal preconditioned RestartedGMResSolver. More...
 
class  Dumux::BlockDiagAMGPreconditioner< M, X, Y, blockLevel >
 A simple ilu0 block diagonal preconditioner. More...
 
class  Dumux::BlockDiagAMGBiCGSTABSolver
 A simple ilu0 block diagonal preconditioned BiCGSTABSolver. More...
 
class  Dumux::LinearSolver
 Base class for linear solvers. More...
 
class  Dumux::Detail::StokesPreconditioner< M, X, Y, l >
 A Stokes preconditioner (saddle-point problem) for the problem \( \begin{pmatrix} A & B \\ C & 0 \end{pmatrix} \begin{pmatrix} u \\ p \end{pmatrix} = \begin{pmatrix} f \\ g \end{pmatrix}, \). More...
 
class  Dumux::StokesSolver< Matrix, Vector, VelocityGG, PressureGG >
 Preconditioned iterative solver for the incompressible Stokes problem. More...
 

Typedefs

template<class ... T>
using Dumux::IstlSolverFactoryBackend = typename Detail::IstlSolverFactoryBackendImplHelper< sizeof...(T)>::template type< T... >
 A linear solver using the dune-istl solver factory to choose the solver and preconditioner at runtime. More...
 
template<class LSTraits , class LATraits >
using Dumux::ILUBiCGSTABIstlSolver = Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::BiCGSTABSolver< typename LATraits::SingleTypeVector >, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory< Dune::SeqILU >, true >
 An ILU preconditioned BiCGSTAB solver using dune-istl. More...
 
template<class LSTraits , class LATraits >
using Dumux::ILURestartedGMResIstlSolver = Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::RestartedGMResSolver< typename LATraits::SingleTypeVector >, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory< Dune::SeqILU >, true >
 An ILU preconditioned GMres solver using dune-istl. More...
 
template<class LSTraits , class LATraits >
using Dumux::SSORBiCGSTABIstlSolver = Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::BiCGSTABSolver< typename LATraits::Vector >, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory< Dune::SeqSSOR > >
 An SSOR-preconditioned BiCGSTAB solver using dune-istl. More...
 
template<class LSTraits , class LATraits >
using Dumux::SSORCGIstlSolver = Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::CGSolver< typename LATraits::Vector >, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory< Dune::SeqSSOR > >
 An SSOR-preconditioned CG solver using dune-istl. More...
 
template<class LSTraits , class LATraits >
using Dumux::AMGBiCGSTABIstlSolver = Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::BiCGSTABSolver< typename LATraits::SingleTypeVector >, Detail::IstlSolvers::IstlAmgPreconditionerFactory, true >
 An AMG preconditioned BiCGSTAB solver using dune-istl. More...
 
template<class LSTraits , class LATraits >
using Dumux::AMGCGIstlSolver = Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::CGSolver< typename LATraits::SingleTypeVector >, Detail::IstlSolvers::IstlAmgPreconditionerFactory, true >
 An AMG preconditioned CG solver using dune-istl. More...
 
template<class LSTraits , class LATraits >
using Dumux::UzawaBiCGSTABIstlSolver = Detail::IstlIterativeLinearSolver< LSTraits, LATraits, Dune::BiCGSTABSolver< typename LATraits::Vector >, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory< Dumux::SeqUzawa > >
 An Uzawa preconditioned BiCGSTAB solver using dune-istl. More...
 
template<class LinearSolverTraits >
using Dumux::ParallelISTLHelper = Detail::ParallelISTLHelperImpl< LinearSolverTraits, LinearSolverTraits::canCommunicate >
 A parallel helper class providing a parallel decomposition of all degrees of freedom. More...
 

Functions

template<class M >
constexpr std::size_t Dumux::Detail::IstlSolvers::preconditionerBlockLevel () noexcept
 Returns the block level for the preconditioner for a given matrix. More...
 
template<class M >
constexpr std::size_t Dumux::preconditionerBlockLevel () noexcept
 Returns the block level for the preconditioner for a given matrix. More...
 
template<class MBlock , class Vector >
void Dumux::symmetrizeConstraints (Dune::BCRSMatrix< MBlock > &A, Vector &b, const Vector &constrainedRows)
 Symmetrize the constrained system Ax=b. More...
 
template<class... MBlock, class Vector >
void Dumux::symmetrizeConstraints (Dune::MultiTypeBlockMatrix< MBlock... > &A, Vector &b, const Vector &constrainedRows)
 Symmetrize the constrained system Ax=b. More...
 

Typedef Documentation

◆ AMGBiCGSTABIstlSolver

template<class LSTraits , class LATraits >
using Dumux::AMGBiCGSTABIstlSolver = typedef Detail::IstlIterativeLinearSolver<LSTraits, LATraits, Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>, Detail::IstlSolvers::IstlAmgPreconditionerFactory, true >

Solver: The BiCGSTAB (stabilized biconjugate gradients method) solver has faster and smoother convergence than the original BiCG. While, it can be applied to nonsymmetric matrices, the preconditioner SSOR assumes symmetry.
See: Van der Vorst, H. A. (1992). "Bi-CGSTAB: A Fast and Smoothly Converging Variant of Bi-CG for the Solution of Nonsymmetric Linear Systems". SIAM J. Sci. and Stat. Comput. 13 (2): 631–644. doi:10.1137/0913035.

Preconditioner: AMG (algebraic multigrid)

◆ AMGCGIstlSolver

template<class LSTraits , class LATraits >
using Dumux::AMGCGIstlSolver = typedef Detail::IstlIterativeLinearSolver<LSTraits, LATraits, Dune::CGSolver<typename LATraits::SingleTypeVector>, Detail::IstlSolvers::IstlAmgPreconditionerFactory, true >

Solver: CG (conjugate gradient) is an iterative method for solving linear systems with a symmetric, positive definite matrix.
See: Helfenstein, R., Koko, J. (2010). "Parallel preconditioned conjugate gradient algorithm on GPU", Journal of Computational and Applied Mathematics, Volume 236, Issue 15, Pages 3584–3590, http://dx.doi.org/10.1016/j.cam.2011.04.025.

Preconditioner: AMG (algebraic multigrid)

◆ ILUBiCGSTABIstlSolver

template<class LSTraits , class LATraits >
using Dumux::ILUBiCGSTABIstlSolver = typedef Detail::IstlIterativeLinearSolver<LSTraits, LATraits, Dune::BiCGSTABSolver<typename LATraits::SingleTypeVector>, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory<Dune::SeqILU>, true >

Solver: The BiCGSTAB (stabilized biconjugate gradients method) solver has faster and smoother convergence than the original BiCG. It can be applied to nonsymmetric matrices.
See: Van der Vorst, H. A. (1992). "Bi-CGSTAB: A Fast and Smoothly Converging Variant of Bi-CG for the Solution of Nonsymmetric Linear Systems". SIAM J. Sci. and Stat. Comput. 13 (2): 631–644. doi:10.1137/0913035.

Preconditioner: ILU(n) incomplete LU factorization. The order n indicates fill-in. It can be damped by the relaxation parameter LinearSolver.PreconditionerRelaxation.
See: Golub, G. H., and Van Loan, C. F. (2012). Matrix computations. JHU Press.

◆ ILURestartedGMResIstlSolver

template<class LSTraits , class LATraits >
using Dumux::ILURestartedGMResIstlSolver = typedef Detail::IstlIterativeLinearSolver<LSTraits, LATraits, Dune::RestartedGMResSolver<typename LATraits::SingleTypeVector>, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory<Dune::SeqILU>, true >

Solver: The GMRes (generalized minimal residual) method is an iterative method for the numerical solution of a nonsymmetric system of linear equations.
See: Saad, Y., Schultz, M. H. (1986). "GMRES: A generalized minimal residual algorithm for solving nonsymmetric linear systems." SIAM J. Sci. and Stat. Comput. 7: 856–869.

Preconditioner: ILU(n) incomplete LU factorization. The order n indicates fill-in. It can be damped by the relaxation parameter LinearSolver.PreconditionerRelaxation.
See: Golub, G. H., and Van Loan, C. F. (2012). Matrix computations. JHU Press.

◆ IstlSolverFactoryBackend

template<class ... T>
using Dumux::IstlSolverFactoryBackend = typedef typename Detail::IstlSolverFactoryBackendImplHelper<sizeof...(T)>::template type<T...>
Note
the solvers are configured via the input file
The template arguments are LinearSolverTraits, LinearAlgebraTraits

◆ ParallelISTLHelper

template<class LinearSolverTraits >
using Dumux::ParallelISTLHelper = typedef Detail::ParallelISTLHelperImpl< LinearSolverTraits, LinearSolverTraits::canCommunicate >

◆ SSORBiCGSTABIstlSolver

template<class LSTraits , class LATraits >
using Dumux::SSORBiCGSTABIstlSolver = typedef Detail::IstlIterativeLinearSolver<LSTraits, LATraits, Dune::BiCGSTABSolver<typename LATraits::Vector>, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory<Dune::SeqSSOR> >

Solver: The BiCGSTAB (stabilized biconjugate gradients method) solver has faster and smoother convergence than the original BiCG. While, it can be applied to nonsymmetric matrices, the preconditioner SSOR assumes symmetry.
See: Van der Vorst, H. A. (1992). "Bi-CGSTAB: A Fast and Smoothly Converging Variant of Bi-CG for the Solution of Nonsymmetric Linear Systems". SIAM J. Sci. and Stat. Comput. 13 (2): 631–644. doi:10.1137/0913035.

Preconditioner: SSOR symmetric successive overrelaxation method. The relaxation is controlled by the parameter LinearSolver.PreconditionerRelaxation. In each preconditioning step, it is applied as often as given by the parameter LinearSolver.PreconditionerIterations.
See: Golub, G. H., and Van Loan, C. F. (2012). Matrix computations. JHU Press.

◆ SSORCGIstlSolver

template<class LSTraits , class LATraits >
using Dumux::SSORCGIstlSolver = typedef Detail::IstlIterativeLinearSolver<LSTraits, LATraits, Dune::CGSolver<typename LATraits::Vector>, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory<Dune::SeqSSOR> >

Solver: CG (conjugate gradient) is an iterative method for solving linear systems with a symmetric, positive definite matrix.
See: Helfenstein, R., Koko, J. (2010). "Parallel preconditioned conjugate gradient algorithm on GPU", Journal of Computational and Applied Mathematics, Volume 236, Issue 15, Pages 3584–3590, http://dx.doi.org/10.1016/j.cam.2011.04.025.

Preconditioner: SSOR symmetric successive overrelaxation method. The relaxation is controlled by the parameter LinearSolver.PreconditionerRelaxation. In each preconditioning step, it is applied as often as given by the parameter LinearSolver.PreconditionerIterations.
See: Golub, G. H., and Van Loan, C. F. (2012). Matrix computations. JHU Press.

◆ UzawaBiCGSTABIstlSolver

template<class LSTraits , class LATraits >
using Dumux::UzawaBiCGSTABIstlSolver = typedef Detail::IstlIterativeLinearSolver<LSTraits, LATraits, Dune::BiCGSTABSolver<typename LATraits::Vector>, Detail::IstlSolvers::IstlDefaultBlockLevelPreconditionerFactory<Dumux::SeqUzawa> >

Solver: The BiCGSTAB (stabilized biconjugate gradients method) solver has faster and smoother convergence than the original BiCG. While, it can be applied to nonsymmetric matrices, the preconditioner SSOR assumes symmetry.
See: Van der Vorst, H. A. (1992). "Bi-CGSTAB: A Fast and Smoothly Converging Variant of Bi-CG for the Solution of Nonsymmetric Linear Systems". SIAM J. Sci. and Stat. Comput. 13 (2): 631–644. doi:10.1137/0913035.

Preconditioner: Uzawa method for saddle point problems

Note
Expects a 2x2 MultiTypeBlockMatrix

Function Documentation

◆ preconditionerBlockLevel() [1/2]

template<class M >
constexpr std::size_t Dumux::Detail::IstlSolvers::preconditionerBlockLevel ( )
constexprnoexcept
Template Parameters
MThe matrix.

◆ preconditionerBlockLevel() [2/2]

template<class M >
constexpr std::size_t Dumux::preconditionerBlockLevel ( )
constexprnoexcept
Template Parameters
MThe matrix.

◆ symmetrizeConstraints() [1/2]

template<class MBlock , class Vector >
void Dumux::symmetrizeConstraints ( Dune::BCRSMatrix< MBlock > &  A,
Vector &  b,
const Vector &  constrainedRows 
)
Note
The function addresses the case where some rows in the linear system have been replaced to fix a degree of freedom to a fixed value. A typical use case is strongly enforced Dirichlet constraints. To symmetrize the constraints in operator A, this function eliminates the corresponding column entries by bringing them to the right-hand side. Particularly, if matrix operator A was symmetric before incorporating constraints, calling this function restores the symmetry of operator A.
the constrainedRows vector is the same shape and type as b and contains 1 for each row that corresponds to constraints (e.g. strongly enforced Dirichlet dofs) and 0 elsewhere
this is a specialization for Dune::BCRSMatrix

◆ symmetrizeConstraints() [2/2]

template<class... MBlock, class Vector >
void Dumux::symmetrizeConstraints ( Dune::MultiTypeBlockMatrix< MBlock... > &  A,
Vector &  b,
const Vector &  constrainedRows 
)
Note
The function addresses the case where some rows in the linear system have been replaced to fix a degree of freedom to a fixed value. A typical use case is strongly enforced Dirichlet constraints. To symmetrize the constraints in operator A, this function eliminates the corresponding column entries by bringing them to the right-hand side. Particularly, if matrix operator A was symmetric before incorporating constraints, calling this function restores the symmetry of operator A.
the constrainedRows vector is the same shape and type as b and contains 1 for each row that corresponds to constraints (e.g. strongly enforced Dirichlet dofs) and 0 elsewhere
this is a specialization for Dune::MultiTypeBlockMatrix