Preconditioner: Multigrid
LFAToolkit.Multigrid — TypeMultigrid(fineoperator, coarseoperator, smoother, prolongation, multigridtype)Multigrid preconditioner for finite element operators
Arguments:
fineoperator::Operator: finite element operator to preconditioncoarseoperator::Union{Operator,Multigrid}: coarse grid representation of finite element operator to preconditionsmoother::AbstractPreconditioner: error relaxation operator, such as Jacobiprolongationbases::AbstractArray{<:AbstractBasis}: element prolongation bases from coarse to fine gridmultigridtype::MultigridType: multigrid type, h-multigrid or p-multigrid
Returns:
- multigrid preconditioner object
LFAToolkit.getnodecoordinatedifferences — Methodgetnodecoordinateddifferences(multigrid)Compute or retrieve the array of differences in coordinates between nodes
Returns:
- array of differences in coordinates between nodes
LFAToolkit.getprolongationmatrix — Functiongetprolongationmatrix(multigrid)Compute or retrieve the prolongation matrix
Returns:
- matrix prolonging from coarse grid to fine grid
Example:
# setup
mesh = Mesh2D(1.0, 1.0);
ctofbasis = TensorH1LagrangeBasis(2, 3, 1, 2; collocatedquadrature = true);
# operators
finediffusion = GalleryOperator("diffusion", 3, 3, mesh);
coarsediffusion = GalleryOperator("diffusion", 2, 3, mesh);
# smoother
jacobi = Jacobi(finediffusion);
# preconditioner
multigrid = PMultigrid(finediffusion, coarsediffusion, jacobi, [ctofbasis]);
# verify
u = ones(ctofbasis.numbernodes);
v = multigrid.prolongationmatrix * u;
@assert v' ≈ [4.0 2.0 4.0 2.0 1.0 2.0 4.0 2.0 4.0] .^ -1
# output
LFAToolkit.computesymbolsprolongation — Methodcomputesymbolsprolongation(multigrid, θ)Compute the symbol matrix for a multigrid prolongation operator
Arguments:
multigrid: Multigrid operator to compute prolongation symbol matrix forθ: Fourier mode frequency array (one frequency per dimension)
Returns:
- symbol matrix for the multigrid prolongation operator