DUNE PDELab (2.8)
Local operator that can be used to create a fully matrix-free Jacobi preconditioner. More...
#include <dune/pdelab/backend/istl/matrixfree/iterativeblockjacobipreconditioner.hh>
Public Types | |
Flags selective assembly | |
enum | |
Whether to do selective assembly on the elements, i.e. whether or not skip_entity() should be called. | |
enum | |
Whether to do selective assembly on the intersections, i.e. whether or not skip_intersection() should be called. | |
Flags for the sparsity pattern | |
enum | |
Whether to assemble the pattern on the elements, i.e. whether or not pattern_volume() should be called. | |
enum | |
Whether to assemble the pattern on the elements after the skeleton has been handled, i.e. whether or not pattern_volume_post_skeleton() should be called. | |
enum | |
Whether to assemble the pattern on the interior intersections, i.e. whether or not pattern_skeleton() should be called. | |
enum | |
Whether to assemble the pattern on the boundary intersections, i.e. whether or not pattern_boundary() should be called. | |
Flags for the non-constant part of the residual and the jacobian | |
enum | |
Whether to call the local operator's alpha_volume(), jacobian_apply_volume() and jacobian_volume(). | |
enum | |
Whether to call the local operator's alpha_volume_post_skeleton(), jacobian_apply_volume_post_skeleton() and jacobian_volume_post_skeleton(). | |
enum | |
Whether to call the local operator's alpha_skeleton(), jacobian_apply_skeleton() and jacobian_skeleton(). | |
enum | |
Whether to call the local operator's alpha_boundary(), jacobian_apply_boundary() and jacobian_boundary(). | |
Flags for the constant part of the residual | |
enum | |
Whether to call the local operator's lambda_volume(). | |
enum | |
Whether to call the local operator's lambda_volume_post_skeleton(). | |
enum | |
Whether to call the local operator's lambda_skeleton(). | |
enum | |
Whether to call the local operator's lambda_boundary(). | |
Special flags | |
enum | |
Whether to visit the skeleton methods from both sides. | |
enum | |
Wheter the local operator describes a linear problem. | |
Public Member Functions | |
IterativeBlockJacobiPreconditionerLocalOperator (const BlockDiagonalLocalOperator &blockDiagonalLocalOperator, const PointDiagonalLocalOperator &pointDiagonalLocalOperator, const GridFunctionSpace &gridFunctionSpace, SolverStatistics< int > &solverStatiscits, BlockSolverOptions solveroptions, const bool verbose=0) | |
Constructor. More... | |
template<typename EG , typename LFSU , typename X , typename LFSV , typename Y > | |
void | alpha_volume (const EG &eg, const LFSU &lfsu, const X &x, const LFSV &lfsv, Y &y) const |
Prepare fully matrix-free preconditioner. More... | |
template<typename EG , typename LFSU , typename Z , typename LFSV , typename Y > | |
void | jacobian_apply_volume (const EG &eg, const LFSU &lfsu, const Z &z, const LFSV &lfsv, Y &y) const |
Apply fully matrix-free preconditioner - linear case. | |
template<typename EG , typename LFSU , typename X , typename Z , typename LFSV , typename Y > | |
void | jacobian_apply_volume (const EG &eg, const LFSU &lfsu, const X &x, const Z &z, const LFSV &lfsv, Y &y) const |
Apply fully matrix-free preconditioner - nonlinear case. More... | |
Detailed Description
class Dune::PDELab::IterativeBlockJacobiPreconditionerLocalOperator< BlockDiagonalLocalOperator, PointDiagonalLocalOperator, GridFunctionSpace, DomainField, IterativeSolver >
Local operator that can be used to create a fully matrix-free Jacobi preconditioner.
Similar to the partial matrix-free class AssembledBlockJacobiPreconditionerLocalOperator this implements a local operator that can be used to implement a matrix-free Jacobi preconditioner. In contrast to the other class this implementation will be fully matrix-free.
A matrix-free Jacobi preconditioner needs to be able to apply the inverse of the block diagonal D to a vector. In the partial matrix-free class mentioned above this was done by assembling the diagonal blocks and inverting this matrix. In order to get a fully matrix-free version we instead use a Krylow solver on the diagonal block. This can once again be done in a matrix-free way but we will need a preconditioner for iterative solver. For this purpose we use another Jacobi preconditioner that operates on a single block. This means we need to assemble the point diagonal of this block and apply the inverse.
For examples see dune-pdelab/dune/pdelab/test/matrixfree/.
- Template Parameters
-
BlockDiagonalLocalOperator A lop for local application of diagonal blocks PointDiagonalLocalOperator A lop for local assembly of point diagonal GridFunctionSpace A grid function space DomainField Domain field IterativeSolver Solver that will be used to 'invert' the diagonal blocks
Constructor & Destructor Documentation
◆ IterativeBlockJacobiPreconditionerLocalOperator()
|
inline |
Constructor.
- Parameters
-
blockDiagonalLocalOperator A local operator that can be used to (locally) apply a diagonal block through the <applyLocalDiagonalBlock> function. You can create such a local operator by wrapping your local operator into a <BlockDiagonalLocalOperatorWrapper> pointDiagonalLocalOperator A local operator that can be used to (locally) assemble the point diagonal of a diagonal block. You can create such a local operator by wrapping your local operator into a <PointDiagonalLocalOperatorWrapper> gridFunctionSpace A grid function space solverStat A class for export solver statistics solveroptions A class for configuring your solver verbose Controls the amount of output
Member Function Documentation
◆ alpha_volume()
|
inline |
Prepare fully matrix-free preconditioner.
This assembles the point diagonal of the diagonal block and stores its inverse. During the apply step this will be used as a preconditioner for the solver on the local block.
◆ jacobian_apply_volume()
|
inline |
Apply fully matrix-free preconditioner - nonlinear case.
Compared to the linear case this needs to set the correct linearization point in the BlockDiagonalOperator
References Dune::PDELab::BlockSolverOptions::_maxiter, Dune::PDELab::BlockSolverOptions::_precondition, Dune::PDELab::BlockSolverOptions::_resreduction, Dune::PDELab::BlockSolverOptions::_verbose, Dune::PDELab::BlockSolverOptions::_weight, Dune::PDELab::SolverStatistics< T >::append(), Dune::IterativeSolver< X, Y >::apply(), Dune::PDELab::LocalVector< T, LFSFlavorTag, W >::data(), Dune::SingleCodimSingleGeomTypeMapper< GV, c >::index(), and Dune::InverseOperatorResult::iterations.
The documentation for this class was generated from the following file:
- dune/pdelab/backend/istl/matrixfree/iterativeblockjacobipreconditioner.hh