DUNE PDELab (git)
Dune::PDELab::ISTLBackend_NOVLP_BCGS_SSORk< GO > Class Template Reference
Nonoverlapping parallel BiCGSTAB solver preconditioned by block SSOR. More...
#include <dune/pdelab/backend/istl/novlpistlsolverbackend.hh>
Public Member Functions | |
ISTLBackend_NOVLP_BCGS_SSORk (const GO &grid_operator, unsigned maxiter_=5000, int steps_=5, int verbose_=1) | |
make a linear solver object More... | |
Vector::ElementType | norm (const Vector &v) const |
Compute global norm of a vector. More... | |
void | apply (M &A, V &z, W &r, typename V::ElementType reduction) |
Solve the given linear system. More... | |
const Dune::PDELab::LinearSolverResult< double > & | result () const |
Return access to result data. | |
Detailed Description
template<class GO>
class Dune::PDELab::ISTLBackend_NOVLP_BCGS_SSORk< GO >
class Dune::PDELab::ISTLBackend_NOVLP_BCGS_SSORk< GO >
Nonoverlapping parallel BiCGSTAB solver preconditioned by block SSOR.
- Template Parameters
-
GO The type of the grid operator used for the spatial discretization (or the fakeGOTraits class for the old grid operator space). It is used to adjust the discretization matrix and extract the trial grid function space.
The solver uses a NonoverlappingBlockPreconditioner with underlying sequential SSOR preconditioner. The crucial step is to add up the matrix entries corresponding to the border vertices on each process. This is achieved by performing a VertexExchanger::sumEntries(Matrix&) before constructing the sequential SSOR.
Constructor & Destructor Documentation
◆ ISTLBackend_NOVLP_BCGS_SSORk()
template<class GO >
|
inlineexplicit |
make a linear solver object
- Parameters
-
[in] gfs_ a grid function space [in] maxiter_ maximum number of iterations to do [in] steps_ number of SSOR steps to apply as inner iteration [in] verbose_ print messages if true
Member Function Documentation
◆ apply()
|
inlineinherited |
Solve the given linear system.
- Parameters
-
[in] A the given matrix [out] z the solution vector to be computed [in] r right hand side [in] reduction to be achieved
◆ norm()
|
inlineinherited |
Compute global norm of a vector.
- Parameters
-
[in] v the given vector
The documentation for this class was generated from the following file:
- dune/pdelab/backend/istl/novlpistlsolverbackend.hh
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Nov 12, 23:30, 2024)