KASKADE 7 development version
|
Classes and functions for multigrid preconditioners and solvers. More...
Classes | |
class | Kaskade::AdditiveMultiGrid< Smoother, Prolongation, CoarsePreconditioner > |
An additive multigrid preconditioner for P1 elements. More... | |
class | Kaskade::BDDC::InterfaceAverages< m, Index > |
A class for computing and providing BDDC coarse space constraints (interface averages). More... | |
class | Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer > |
A class representing small KKT systems on a single subdomain. More... | |
class | Kaskade::BDDC::SharedMemoryDomain< Subdomain > |
A class that orchestrates the exchange of data between BDDC subdomains. More... | |
class | Kaskade::BDDC::KKTSolver< Domain > |
A class for solving a block-diagonal-structured KKT system. More... | |
class | Kaskade::BDDC::BDDCSolver< Subdomain > |
Balanced Domain Decomposition with Constraints preconditioner. More... | |
class | Kaskade::BDDC::SpaceTransferBase< m, Scalar > |
Base class for BDDC SpaceTransfer implementations. More... | |
class | Kaskade::BDDC::SpaceTransfer< m, Scalar, TransmissionScalar > |
A trivial BDDC space transfer. More... | |
class | Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation > |
A general multiplicative multigrid preconditioner. More... | |
class | Kaskade::MGProlongation |
A prolongation operator for P1 finite elements from a coarser grid level to the next finer level. More... | |
class | Kaskade::MultiGridStack< Prolongation, Entry, Index > |
Class for multigrid stacks. More... | |
Typedefs | |
using | Kaskade::BDDC::Interfaces = std::vector< std::pair< int, std::vector< int > > > |
A data structure describing the interfaces of a subdomain. More... | |
Enumerations | |
enum | Kaskade::BDDC::InterfaceType { Kaskade::BDDC::CORNER = 1 , Kaskade::BDDC::EDGE = 2 , Kaskade::BDDC::FACE = 4 } |
Discrimination of different interface types. More... | |
Functions | |
template<class Index > | |
std::vector< std::vector< Index > > | Kaskade::BDDC::algebraicDomainDecomposition (NumaCRSPattern< Index > const &A, int n) |
Creates a decomposition of indices into overlapping subsets that can be used for constructing nonoverlapping domain decomposition solvers. More... | |
template<class Scalar , int m, class Index > | |
std::tuple< std::vector< NumaBCRSMatrix< Dune::FieldMatrix< Scalar, m, m >, Index > >, std::vector< Dune::BlockVector< Dune::FieldVector< Scalar, m > > >, std::vector< std::vector< Index > >, std::vector< std::vector< LocalDof > > > | Kaskade::BDDC::matrixDomainDecommposition (NumaBCRSMatrix< Dune::FieldMatrix< Scalar, m, m >, Index > const &A, Dune::BlockVector< Dune::FieldVector< Scalar, m > > const &f, int n, Scalar tolerance=1000 *std::numeric_limits< Scalar >::epsilon(), int maxIter=100) |
template<class Entry , class Index , class Prolongation , class MakeSmoother , class MakeCoarsePreconditioner > | |
auto | Kaskade::makeAdditiveMultiGrid (NumaBCRSMatrix< Entry, Index > const &A, std::vector< Prolongation > const &prolongations, MakeSmoother const &makeSmoother, MakeCoarsePreconditioner const &makeCoarsePreconditioner, bool onlyLowerTriangle=false) |
Convenience function creating additive multigrid preconditioners for P1 elements. More... | |
template<class Entry , class Index , class GridMan > | |
auto | Kaskade::makeBPX (NumaBCRSMatrix< Entry, Index > const &A, GridMan const &gridman, bool onlyLowerTriangle=false) |
Convenience function creating a BPX preconditioner for P1 elements. More... | |
template<class Entry , class Index , class FineSpace , class MakeSmoother , class MakeCoarseSmoother > | |
auto | Kaskade::makeAdditivePMultiGrid (NumaBCRSMatrix< Entry, Index > A, FineSpace const &space, H1Space< typename FineSpace::Grid > const &p1Space, MakeSmoother const &makeSmoother, MakeCoarseSmoother const &makeCoarseSmoother, int coarseLevel=0, bool onlyLowerTriangle=false) |
Convenience function creating additive multigrid preconditioner of V-cycle type for higher order elements. More... | |
template<class StorageTag , class Entry , class Index , class FineSpace > | |
auto | Kaskade::makePBPX (NumaBCRSMatrix< Entry, Index > A, FineSpace const &space, H1Space< typename FineSpace::Grid > const &p1Space, StorageTag storageTag, int coarseLevel=0, bool onlyLowerTriangle=false) |
Convenience function creating additive multigrid preconditioner of V-cycle type with domain decomposition and Jacobi smoother for higher order elements. More... | |
template<class Entry , class Index , class Space > | |
auto | Kaskade::makeMultigrid (NumaBCRSMatrix< Entry, Index > const &A, Space const &space, bool onlyLowerTriangle=false, DirectType directType=DirectType::MUMPS) |
Convenience function for creating a (multiplicative) multigrid preconditioner for elliptic equations. More... | |
template<typename Matrix > | |
auto | Kaskade::makeDirectPreconditioner (Matrix &&A, DirectType directType=DirectType::MUMPS) |
Creates a direct solver for the given matrix. More... | |
template<int domainBlockSize, int rangeBlockSize, class Index > | |
auto | Kaskade::makeDirectPreconditioner (NumaBCRSMatrix< Dune::FieldMatrix< float, rangeBlockSize, domainBlockSize >, Index > &&A, DirectType directType=DirectType::MUMPS) |
Creates a direct solver for the given matrix. This specialization is for matrices containing floats (single precision) as scalar entries. Since the Kaskade interface for direct solvers only works for doubles, the matrix is copied into a new matrix with double as scalar entries. Since the coarse grid matrix often is much smaller, this hopefully means neglectable overhead. More... | |
template<class SparseIndex = size_t, class CoarseSpace , class FineSpace > | |
NumaBCRSMatrix< Dune::FieldMatrix< typename FineSpace::Scalar, 1, 1 >, SparseIndex > | Kaskade::prolongation (CoarseSpace const &coarseSpace, FineSpace const &fineSpace) |
Computes an interpolation-based prolongation matrix from a (supposedly) coarser space to a finer space. More... | |
template<class SparseIndex , class Mapper > | |
std::vector< NumaBCRSMatrix< Dune::FieldMatrix< typename Mapper::Scalar, 1, 1 >, SparseIndex > > | Kaskade::prolongationStack (FEFunctionSpace< Mapper > const &space) |
Computes a stack of prolongation matrices for higher order finite element spaces. More... | |
template<class Entry , class Index > | |
auto | Kaskade::conjugation (MGProlongation const &p, NumaBCRSMatrix< Entry, Index > const &a, bool onlyLowerTriangle=false) |
Creates a Galerkin projected Matrix \( P^T A P \) from a prolongation \( P \) and a symmetric matrix \( A \). More... | |
template<class GridMan > | |
std::vector< MGProlongation > | Kaskade::prolongationStack (GridMan const &gridman, int coarseLevel=0, size_t minNodes=0) |
Computes a sequence of prolongation matrices for P1 finite elements in hierarchical grids. More... | |
template<class GridMan > | |
std::vector< NumaBCRSMatrix< Dune::FieldMatrix< double, 1, 1 > > > | Kaskade::deformedProlongationStack (GridMan const &gridman, int coarseLevel=0, size_t minNodes=0, bool interpolateAtShiftedPosition=true) |
Computes a sequence of prolongation matrices for P1 finite elements in hierarchical grids. More... | |
template<class GridMan , class Entry , class Index > | |
MultiGridStack< MGProlongation, Entry, Index > | Kaskade::makeGeometricMultiGridStack (GridMan const &gridManager, NumaBCRSMatrix< Entry, Index > &&A, size_t minNodes=10000, bool onlyLowerTriangle=false) |
convenience routine for creating multigrid stacks based on geometric coarsening for P1 elements More... | |
template<class Entry , class Index > | |
MultiGridStack< MGProlongation, Entry, Index > | Kaskade::makeAlgebraicMultigridStack (NumaBCRSMatrix< Entry, Index > &&A, Index n=0, bool onlyLowerTriangle=false) |
Creates stack of prolongations and projected Galerkin matrices. More... | |
template<typename FineSpace > | |
auto | Kaskade::makeDeepProlongationStack (FineSpace const &space, int coarseLevel=0, int minNodes=0) |
Convenience routine for creating a deep prolongation stack from coarse to fine grid and on from P1 to higher order. More... | |
template<typename FineSpace , typename Matrix > | |
auto | Kaskade::makePMultiGridStack (FineSpace const &space, Matrix &&A, bool onlyLowerTriangle) |
convenience routine for creating multigrid stacks based on coarsening by reducing the ansatz order from P to P1 More... | |
template<typename FineSpace , typename CoarseSpace , typename Matrix > | |
auto | Kaskade::makePMultiGridStack (FineSpace const &fineSpace, Matrix &&fA, CoarseSpace const &coarseSpace) |
convenience routine for creating multigrid stacks between two spaces. More... | |
template<typename FineSpace , typename CoarseSpace , typename Matrix > | |
auto | Kaskade::makePMultiGridStack (FineSpace const &fineSpace, Matrix &&fA, CoarseSpace const &coarseSpace, Matrix &&cA) |
convenience routine for creating multigrid stacks between two spaces. More... | |
template<class Entry , class Index , class Prolongation , class MakeSmoother , class CoarsePreconditioner > | |
auto | makeMultiplicativeMultiGrid (MultiGridStack< Prolongation, Entry, Index > &&mgStack, MakeSmoother const &makeSmoother, std::unique_ptr< CoarsePreconditioner > &&coarsePreconditioner, int nPre=3, int nPost=3) |
Convenience function creating multiplicative multigrid preconditioner of V-cycle type for P1 elements. More... | |
template<class Entry , class Index , class GridMan > | |
auto | makeJacobiMultiGrid (NumaBCRSMatrix< Entry, Index > const &A, GridMan const &gridman, int nPre=3, int nPost=3, bool onlyLowerTriangle=false, DirectType directType=DirectType::MUMPS) |
Convenience function creating multiplicative Jacobi multigrid for P1 elements. More... | |
template<class Entry , class Index , class Space , class MakeSmoother > | |
auto | makeMultiplicativePMultiGrid (NumaBCRSMatrix< Entry, Index > &&A, Space const &space, MakeSmoother const &makeSmoother, int nPre=3, int nPost=3, bool onlyLowerTriangle=false, DirectType directType=DirectType::MUMPS) |
Convenience function creating multiplicative multigrid preconditioner of V-cycle type for higher order elements. More... | |
template<class Entry , class Index , class Space , class CoarseSpace > | |
auto | makeBlockJacobiPMultiGrid (NumaBCRSMatrix< Entry, Index > A, Space const &space, NumaBCRSMatrix< Entry, Index > coarseA, CoarseSpace const &coarseSpace, int nPre=3, int nPost=3, DirectType directType=DirectType::MUMPS) |
Convenience function creating multiplicative multigrid preconditioner of V-cycle type for higher order elements. More... | |
template<class Entry , class Index , class Space > | |
auto | makeBlockJacobiPMultiGrid (NumaBCRSMatrix< Entry, Index > A, Space const &space, int nPre=3, int nPost=3, bool onlyLowerTriangle=false, DirectType directType=DirectType::MUMPS) |
Convenience function creating multiplicative multigrid preconditioner of V-cycle type for higher order elements. More... | |
template<class Entry , class Index , class FineSpace > | |
auto | makeJacobiPMultiGrid (NumaBCRSMatrix< Entry, Index > A, FineSpace const &space, int nPre=3, int nPost=3, bool onlyLowerTriangle=false, DirectType directType=DirectType::MUMPS) |
Convenience function creating multiplicative multigrid preconditioner of V-cycle type with Jacobi smoother for higher order elements. More... | |
template<class Prolongation , class Entry , class Index > | |
MultiGridStack< Prolongation, Entry, Index > | makeMultiGridStack (std::vector< Prolongation > &&ps, NumaBCRSMatrix< Entry, Index > &&A, bool onlyLowerTriangle) |
Convenience routine for creating multigrid stacks. More... | |
Classes and functions for multigrid preconditioners and solvers.
This includes general multiplicative V-cycles and additive multigrid methods, as well as several convenience functions for defining concrete preconditioners by specifying particular ingredients (smoothers and coarse grid preconditioners, in particular). Multigrid methods for arbitrary finite element order on simplicial grids are available.
If unsure which multigrid method to create, try the default construction makeMultigrid, which tries to come up with reasonable ingredients.
using Kaskade::BDDC::Interfaces = typedef std::vector<std::pair<int,std::vector<int> >> |
A data structure describing the interfaces of a subdomain.
It consists of a sequence of pairs \( (i,idx) \) of adjacent subdomain number \( i \) and the (local to the current subdomain) degrees of freedom (sorted in increasing order) that are shared with subdomain \( i \). The sequence shall be sorted by increasing subdomain number.
Definition at line 34 of file bddcSpaceTransfer.hh.
std::vector< std::vector< Index > > Kaskade::BDDC::algebraicDomainDecomposition | ( | NumaCRSPattern< Index > const & | A, |
int | n | ||
) |
Creates a decomposition of indices into overlapping subsets that can be used for constructing nonoverlapping domain decomposition solvers.
Index | the integral number type to be used for indices |
A | the symmetric sparse matrix pattern |
n | the number of subdomains to create |
Referenced by Kaskade::BDDC::matrixDomainDecommposition().
auto Kaskade::conjugation | ( | MGProlongation const & | p, |
NumaBCRSMatrix< Entry, Index > const & | a, | ||
bool | onlyLowerTriangle = false |
||
) |
Creates a Galerkin projected Matrix \( P^T A P \) from a prolongation \( P \) and a symmetric matrix \( A \).
onlyLowerTriangle | if true, A contains onl the lower triangular part of the symmetric matrix. |
Definition at line 421 of file prolongation.hh.
std::vector< NumaBCRSMatrix< Dune::FieldMatrix< double, 1, 1 > > > Kaskade::deformedProlongationStack | ( | GridMan const & | gridman, |
int | coarseLevel = 0 , |
||
size_t | minNodes = 0 , |
||
bool | interpolateAtShiftedPosition = true |
||
) |
Computes a sequence of prolongation matrices for P1 finite elements in hierarchical grids.
We assume that each vertex not contained in the coarse grid has been created by bisecting an edge, and possibly shifted a moderate distance.
If a vertex is shifted, interpolating their direct parent nodes leads to nonlinear interpolation, i.e. coarse grid functions which are not piecwise linear on the (deformed) coarse grid cells. This can lead to high energy of the coarse grid cells. We therefore interpolate the vertex value not only from its direct parents, but from all corners of a father cell.
Grid | the Dune grid type on which the P1 space is defined. The grid has to be a simplicial grid. |
grid | the grid itself |
coarseLevel | the grid level to use as coarsest level |
minNodes | minimum number of nodes to keep as coarse grid |
interpolateAtShiftedPosition | whether or not potential vertex shifts shall be respected or not |
Definition at line 644 of file prolongation.hh.
Referenced by Kaskade::makeDeepProlongationStack().
auto Kaskade::makeAdditiveMultiGrid | ( | NumaBCRSMatrix< Entry, Index > const & | A, |
std::vector< Prolongation > const & | prolongations, | ||
MakeSmoother const & | makeSmoother, | ||
MakeCoarsePreconditioner const & | makeCoarsePreconditioner, | ||
bool | onlyLowerTriangle = false |
||
) |
Convenience function creating additive multigrid preconditioners for P1 elements.
A | the sparse Galerkin matrix for P1 finite elements on the leaf view of the given simplicial grid |
grid | a simplicial grid |
makeSmoother | a callable object taking a (projected) Galerkin matrix of type NumaBCRSMatrix<Entry,Index> and a level, returning a preconditioner |
onlyLowerTriangle | if true, A is assumed to be symmetric and only its lower triangular part is accessed |
Definition at line 159 of file additiveMultigrid.hh.
Referenced by Kaskade::makeAdditivePMultiGrid(), and Kaskade::makeBPX().
auto Kaskade::makeAdditivePMultiGrid | ( | NumaBCRSMatrix< Entry, Index > | A, |
FineSpace const & | space, | ||
H1Space< typename FineSpace::Grid > const & | p1Space, | ||
MakeSmoother const & | makeSmoother, | ||
MakeCoarseSmoother const & | makeCoarseSmoother, | ||
int | coarseLevel = 0 , |
||
bool | onlyLowerTriangle = false |
||
) |
Convenience function creating additive multigrid preconditioner of V-cycle type for higher order elements.
This is realized as an additive h-multigrid nested within a two-grid scheme. The outer two-level algorithm uses the smoother obtained from
Note that a pure Jacobi preconditioner becomes quickly inefficient when the polynomial ansatz order increases. From degree 3 on, using an overlapping domain decomposition smoother such as PatchDomainDecompositionPreconditioner is probably more efficient.
A | the sparse Galerkin matrix for higher order finite elements on the leaf view of the given simplicial grid |
space | a higher order finite element space |
p1Space | a linear finite element space |
makeSmoother | a callable object taking a (projected) Galerkin matrix of type NumaBCRSMatrix<Entry,Index> and a level, returning a preconditioner |
makeCoarseSmoother | a callable object taking a (projected) Galerkin matrix of type NumaBCRSMatrix<Entry,Index> and a level, returning a preconditioner |
coarseLevel | the grid level down to which the additive h-multigrid shall reach |
onlyLowerTriangle | if true, A is assumed to be symmetric and only its lower triangular part is accessed |
Definition at line 216 of file additiveMultigrid.hh.
Referenced by Kaskade::makePBPX().
MultiGridStack< MGProlongation, Entry, Index > Kaskade::makeAlgebraicMultigridStack | ( | NumaBCRSMatrix< Entry, Index > && | A, |
Index | n = 0 , |
||
bool | onlyLowerTriangle = false |
||
) |
Creates stack of prolongations and projected Galerkin matrices.
A | the symmetric sparse matrix |
n | stop the coarsening if the number of rows/cols drops below this number |
onlyLowerTriangle | if true, only the lower triangular part of symmetric A is accessed |
The prolongation/restriction used is particularly simple: Stepping through the nodes, every not yet classified node is classified as coarse node, and its direct neighbours as fine nodes. Thus every fine node has at least one coarse neighbour, and the coarse nodes are not adjacent.
For each fine node, the at most two neighbouring coarse nodes with largest connection strength are selected as parents. From these, the fine node is interpolated.
|
related |
Convenience function creating multiplicative multigrid preconditioner of V-cycle type for higher order elements.
This is realized as a multiplicative h-multigrid nested within a two-grid scheme. The outer two-level algorithm uses an overlapping domain decomposition smoother on the patches around grid vertices. Classical h-multigrid with point Jacobi smoother is used as a coarse level preconditioner/solver.
The default smoother step sizes are 0.5 for the domain decomposition smoother and 0.8 for the Jacobi smoother. In a preliminary test (Poisson equation on the unit square, 2016-01) these values turned out to be quite good. But for 3D 0.4 seems to be better in elastomechanics with flat geometry.
A | the sparse Galerkin matrix for higher order finite elements on the leaf view of the given simplicial grid |
space | a higher order finite element space |
nPre | the number of pre-smoothings |
nPost | the number of post-smoothings |
onlyLowerTriangle | if true, A is assumed to be symmetric and only its lower triangular part is accessed |
Definition at line 581 of file multiplicativeMultigrid.hh.
|
related |
Convenience function creating multiplicative multigrid preconditioner of V-cycle type for higher order elements.
This is realized as a multiplicative h-multigrid nested within a two-grid scheme. The outer two-level algorithm uses an overlapping domain decomposition smoother on the patches around grid vertices. Classical h-multigrid with point Jacobi smoother is used as a coarse level preconditioner/solver.
The default smoother step sizes are 0.5 for the domain decomposition smoother and 0.8 for the Jacobi smoother. In a preliminary test (Poisson equation on the unit square, 2016-01) these values turned out to be quite good. But for 3D 0.4 seems to be better in elastomechanics with flat geometry.
A | the sparse Galerkin matrix for higher order finite elements on the leaf view of the given simplicial grid |
space | a higher order finite element space |
nPre | the number of pre-smoothings |
nPost | the number of post-smoothings |
onlyLowerTriangle | if true, A is assumed to be symmetric and only its lower triangular part is accessed |
Definition at line 537 of file multiplicativeMultigrid.hh.
auto Kaskade::makeBPX | ( | NumaBCRSMatrix< Entry, Index > const & | A, |
GridMan const & | gridman, | ||
bool | onlyLowerTriangle = false |
||
) |
Convenience function creating a BPX preconditioner for P1 elements.
A | the sparse Galerkin matrix for P1 finite elements on the leaf view of the given simplicial grid |
grid | a simplicial grid |
onlyLowerTriangle | if true, A is assumed to be symmetric and only its lower triangular part is accessed |
Definition at line 180 of file additiveMultigrid.hh.
auto Kaskade::makeDeepProlongationStack | ( | FineSpace const & | space, |
int | coarseLevel = 0 , |
||
int | minNodes = 0 |
||
) |
Convenience routine for creating a deep prolongation stack from coarse to fine grid and on from P1 to higher order.
This relies on a geometric grid hierarchy for the P1 substack.
The multigrid recursion towards coarser levels is terminated as soon as the given coarse level is reached or if the given number of nodes is no longer reached.
Definition at line 921 of file prolongation.hh.
Referenced by Kaskade::makeJointPstack().
auto Kaskade::makeDirectPreconditioner | ( | Matrix && | A, |
DirectType | directType = DirectType::MUMPS |
||
) |
Creates a direct solver for the given matrix.
Definition at line 383 of file multiplicativeMultigrid.hh.
Referenced by Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::makeJacobiMultiGrid().
auto Kaskade::makeDirectPreconditioner | ( | NumaBCRSMatrix< Dune::FieldMatrix< float, rangeBlockSize, domainBlockSize >, Index > && | A, |
DirectType | directType = DirectType::MUMPS |
||
) |
Creates a direct solver for the given matrix. This specialization is for matrices containing floats (single precision) as scalar entries. Since the Kaskade interface for direct solvers only works for doubles, the matrix is copied into a new matrix with double as scalar entries. Since the coarse grid matrix often is much smaller, this hopefully means neglectable overhead.
Definition at line 436 of file multiplicativeMultigrid.hh.
MultiGridStack< MGProlongation, Entry, Index > Kaskade::makeGeometricMultiGridStack | ( | GridMan const & | gridManager, |
NumaBCRSMatrix< Entry, Index > && | A, | ||
size_t | minNodes = 10000 , |
||
bool | onlyLowerTriangle = false |
||
) |
convenience routine for creating multigrid stacks based on geometric coarsening for P1 elements
Definition at line 882 of file prolongation.hh.
Referenced by Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::makeJacobiMultiGrid().
|
related |
Convenience function creating multiplicative Jacobi multigrid for P1 elements.
This creates a V-cycle multigrid preconditioner with Jacobi smoother. A direct solver is used for grid level 0.
A | the sparse Galerkin matrix for P1 finite elements on the leaf view of the given simplicial grid |
gridman | a grid manager base of a simplicial grid |
nPre | the number of pre-smoothings |
nPost | the number of post-smoothings |
onlyLowerTriangle | if true, A is assumed to be symmetric and only its lower triangular part is accessed |
Definition at line 460 of file multiplicativeMultigrid.hh.
|
related |
Convenience function creating multiplicative multigrid preconditioner of V-cycle type with Jacobi smoother for higher order elements.
A | the sparse Galerkin matrix for higher order finite elements on the leaf view of the given simplicial grid |
space | a higher order finite element space |
nPre | the number of pre-smoothings |
nPost | the number of post-smoothings |
onlyLowerTriangle | if true, A is assumed to be symmetric and only its lower triangular part is accessed |
Definition at line 603 of file multiplicativeMultigrid.hh.
auto Kaskade::makeMultigrid | ( | NumaBCRSMatrix< Entry, Index > const & | A, |
Space const & | space, | ||
bool | onlyLowerTriangle = false , |
||
DirectType | directType = DirectType::MUMPS |
||
) |
Convenience function for creating a (multiplicative) multigrid preconditioner for elliptic equations.
This function chooses sensible default values for the structure of the preconditioner depending in particular on the finite element ansatz order and the grid refinement depth. As the resulting preconditioners have different types depending on runtime information, we return a unique pointer to SymmetricPreconditioner. Use this if you don't have more specific ideas on how you would like the multigrid to work.
The function assumes the space consists of globally continuous finite element functions on simplicial grids with scalar shape functions.
Definition at line 40 of file multigrid.hh.
|
related |
Convenience routine for creating multigrid stacks.
Given a stack of prolongations and the top level Galerkin matrix, this creates the complete stack including all projected Galerkin matrices.
ps | a vector of prolongation matrices |
A | the top level (finest) Galerkin matrix to be projected down |
onlyLowerTriangle | if true, only the lower triangular part of A will be referenced |
Definition at line 871 of file prolongation.hh.
|
related |
Convenience function creating multiplicative multigrid preconditioner of V-cycle type for P1 elements.
mgStack | a stack of prolongations and matching projected Galerkin matrices |
makeSmoother | a callable object taking a (projected) Galerkin matrix of type NumaBCRSMatrix<Entry,Index> and a level, returning a preconditioner |
coarsePreconditioner | a symmetric preconditioner to be used for "solving" on the coarsest level |
nPre | the number of pre-smoothings |
nPost | the number of post-smoothings |
Definition at line 320 of file multiplicativeMultigrid.hh.
|
related |
Convenience function creating multiplicative multigrid preconditioner of V-cycle type for higher order elements.
A | the sparse Galerkin matrix for higher order finite elements on the leaf view of the given simplicial grid |
space | a higher order finite element space |
p1Space | a linear finite element space |
makeSmoother | a callable object taking a (projected) Galerkin matrix of type NumaBCRSMatrix<Entry,Index> and a level, returning a preconditioner |
nPre | the number of pre-smoothings |
nPost | the number of post-smoothings |
onlyLowerTriangle | if true, A is assumed to be symmetric and only its lower triangular part is accessed |
Definition at line 496 of file multiplicativeMultigrid.hh.
auto Kaskade::makePBPX | ( | NumaBCRSMatrix< Entry, Index > | A, |
FineSpace const & | space, | ||
H1Space< typename FineSpace::Grid > const & | p1Space, | ||
StorageTag | storageTag, | ||
int | coarseLevel = 0 , |
||
bool | onlyLowerTriangle = false |
||
) |
Convenience function creating additive multigrid preconditioner of V-cycle type with domain decomposition and Jacobi smoother for higher order elements.
StorageTag | a tag type for selecting the storage of inverse local matrices in the overlapping Schwarz smoother |
A | the sparse Galerkin matrix for higher order finite elements on the leaf view of the given simplicial grid |
space | a higher order finite element space |
p1Space | a linear finite element space |
coarseLevel | the grid level down to which the additive h-multigrid shall reach |
onlyLowerTriangle | if true, A is assumed to be symmetric and only its lower triangular part is accessed |
The space is decomposed in (i) overlapping patches, (ii) a P1 mesh hierarchy, and (iii) a coarse space. For (i), direct solvers are used, resulting in an overlapping additive Schwarz smoother. For (ii), standard point Jacobi with damping factor 0.6 is used. (iii) is again treated with a direct solver.
Definition at line 254 of file additiveMultigrid.hh.
auto Kaskade::makePMultiGridStack | ( | FineSpace const & | fineSpace, |
Matrix && | fA, | ||
CoarseSpace const & | coarseSpace | ||
) |
convenience routine for creating multigrid stacks between two spaces.
Definition at line 989 of file prolongation.hh.
auto Kaskade::makePMultiGridStack | ( | FineSpace const & | fineSpace, |
Matrix && | fA, | ||
CoarseSpace const & | coarseSpace, | ||
Matrix && | cA | ||
) |
convenience routine for creating multigrid stacks between two spaces.
An approximation of the projected Galerkin matrix is provided. Often, this can be assembled much more efficiently than projected.
Definition at line 1002 of file prolongation.hh.
Referenced by Kaskade::makePMultiGridStack().
auto Kaskade::makePMultiGridStack | ( | FineSpace const & | space, |
Matrix && | A, | ||
bool | onlyLowerTriangle | ||
) |
convenience routine for creating multigrid stacks based on coarsening by reducing the ansatz order from P to P1
Definition at line 946 of file prolongation.hh.
Referenced by Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::makeBlockJacobiPMultiGrid(), and Kaskade::MultiplicativeMultiGrid< Entry, Index, Smoother, Prolongation >::makeMultiplicativePMultiGrid().
std::tuple< std::vector< NumaBCRSMatrix< Dune::FieldMatrix< Scalar, m, m >, Index > >, std::vector< Dune::BlockVector< Dune::FieldVector< Scalar, m > > >, std::vector< std::vector< Index > >, std::vector< std::vector< LocalDof > > > Kaskade::BDDC::matrixDomainDecommposition | ( | NumaBCRSMatrix< Dune::FieldMatrix< Scalar, m, m >, Index > const & | A, |
Dune::BlockVector< Dune::FieldVector< Scalar, m > > const & | f, | ||
int | n, | ||
Scalar | tolerance = 1000*std::numeric_limits<Scalar>::epsilon() , |
||
int | maxIter = 100 |
||
) |
NumaBCRSMatrix< Dune::FieldMatrix< typename FineSpace::Scalar, 1, 1 >, SparseIndex > Kaskade::prolongation | ( | CoarseSpace const & | coarseSpace, |
FineSpace const & | fineSpace | ||
) |
Computes an interpolation-based prolongation matrix from a (supposedly) coarser space to a finer space.
CoarseSpace | a FEFunctionSpace type |
FineSpace | a FEFunctionSpace type |
coarseSpace | the (supposedly) coarser space (domain) |
fineSpace | the (supposedly) finer space (image) |
Both function spaces have to be defined on the very same grid view and to be geometrically nested.
Definition at line 43 of file prolongation.hh.
Referenced by Kaskade::TaylorHoodPreconditioner< X >::apply(), Kaskade::localProlongationMatrix(), and Kaskade::TransferData< Space, CoarseningPolicy >::transferMatrix().
std::vector< NumaBCRSMatrix< Dune::FieldMatrix< typename Mapper::Scalar, 1, 1 >, SparseIndex > > Kaskade::prolongationStack | ( | FEFunctionSpace< Mapper > const & | space | ) |
Computes a stack of prolongation matrices for higher order finite element spaces.
This computes a stack of prolongation matrices for a scale of finite element spaces with increasing polynomial ansatz order. From level to level, the ansatz order is doubled (maybe +1 for odd degrees).
Mapper | a finite element local to global mapper with scalar shape functions. |
Definition at line 131 of file prolongation.hh.
Referenced by Kaskade::makeAdditivePMultiGrid(), Kaskade::makeBPX(), Kaskade::makeGeometricMultiGridStack(), and Kaskade::prolongationStack().
std::vector< MGProlongation > Kaskade::prolongationStack | ( | GridMan const & | gridman, |
int | coarseLevel = 0 , |
||
size_t | minNodes = 0 |
||
) |
Computes a sequence of prolongation matrices for P1 finite elements in hierarchical grids.
We assume that each vertex not contained in the coarse grid has been created by bisecting an edge, such that there are exactly two "parent vertices".
Grid | the Dune grid type on which the P1 space is defined. The grid has to be a simplicial grid. |
grid | the grid itself |
minNodes | minimum number of nodes to keep as coarse grid |
Definition at line 588 of file prolongation.hh.