Preconditioner: BDDC

LFAToolkit supports lumped and Dirichlet BDDC preconditioners.

LFAToolkit.BDDCInjectionType.BDDCInjectTypeType

BDDC injection types

Types:

  • scaled: scaled injection
  • harmonic: discrete harmonic extension

Example

LFAToolkit.BDDCInjectionType.BDDCInjectType

# output
Enum LFAToolkit.BDDCInjectionType.BDDCInjectType:
scaled = 0
harmonic = 1
source

Lumped BDDC

Examples

This is an example of a simple lumped BDDC preconditioner.

# ------------------------------------------------------------------------------
# lumped BDDC example
# ------------------------------------------------------------------------------

using LFAToolkit
using LinearAlgebra

# setup
mesh = Mesh2D(1.0, 1.0)
p = 3

# diffusion operator
diffusion = GalleryOperator("diffusion", p + 1, p + 1, mesh)

# lumped BDDC preconditioner
bddc = LumpedBDDC(diffusion)

# compute operator symbols
A = computesymbols(bddc, [0.2], [π, π])
eigenvalues = real(eigvals(A))

# ------------------------------------------------------------------------------

This is an example of a simple lumped BDDC preconditioner on a macro-element patch.

# ------------------------------------------------------------------------------
# lumped BDDC on macro elements example
# ------------------------------------------------------------------------------

using LFAToolkit
using LinearAlgebra

# setup
mesh = Mesh2D(1.0, 1.0)
p = 1
numberelements1d = 4

# operator
diffusion = GalleryMacroElementOperator("diffusion", p + 1, p + 2, numberelements1d, mesh)

# BDDC smoother
bddc = LumpedBDDC(diffusion)

# compute operator symbols
A = computesymbols(bddc, [0.2], [π, π])
eigenvalues = real(eigvals(A))

# ------------------------------------------------------------------------------

Documentation

LFAToolkit.LumpedBDDCFunction
LumpedBDDC(operator)

Lumped BDDC preconditioner for finite element operators

Arguments:

  • operator::Operator: finite element operator to precondition

Returns:

  • lumped BDDC preconditioner object

Example:

# setup
mesh = Mesh2D(1.0, 1.0);

# operator
finediffusion = GalleryOperator("diffusion", 5, 5, mesh);

# preconditioner
bddc = LumpedBDDC(finediffusion);

# verify
println(bddc)
println(bddc.operator)

# output

lumped BDDC preconditioner
finite element operator:
2d mesh:
    dx: 1.0
    dy: 1.0

2 inputs:
operator field:
  tensor product basis:
    numbernodes1d: 5
    numberquadraturepoints1d: 5
    numbercomponents: 1
    dimension: 2
  evaluation mode:
    gradient
operator field:
  tensor product basis:
    numbernodes1d: 5
    numberquadraturepoints1d: 5
    numbercomponents: 1
    dimension: 2
  evaluation mode:
    quadratureweights

1 output:
operator field:
  tensor product basis:
    numbernodes1d: 5
    numberquadraturepoints1d: 5
    numbercomponents: 1
    dimension: 2
  evaluation mode:
    gradient
source
LFAToolkit.computesymbolsMethod
computesymbols(bddc, ω, θ)

Compute the symbol matrix for a BDDC preconditioned operator

Arguments:

  • bddc::BDDC: BDDC preconditioner to compute symbol matrix for
  • ω::Array{Real}: BDDC smoothing parameter
  • θ::Array{Real}: Fourier mode frequency array (one frequency per dimension)

Returns:

  • symbol matrix for the BDDC preconditioned operator

Lumped BDDC Example:

using LinearAlgebra;

for dimension = 2:3
    # setup
    mesh = []
    if dimension == 2
        mesh = Mesh2D(1.0, 1.0)
    elseif dimension == 3
        mesh = Mesh3D(1.0, 1.0, 1.0)
    end
    diffusion = GalleryOperator("diffusion", 3, 3, mesh)
    bddc = LumpedBDDC(diffusion)

    # compute symbols
    A = computesymbols(bddc, [0.2], π * ones(dimension))

    # verify
    eigenvalues = real(eigvals(A))
    if dimension == 2
        @assert minimum(eigenvalues) ≈ 0.43999999999999995
        @assert maximum(eigenvalues) ≈ 0.8
    elseif dimension == 3
        @assert minimum(eigenvalues) ≈ -0.6319999999999972
        @assert maximum(eigenvalues) ≈ 0.8
    end
end

# output

Dirichlet BDDC Example:

using LinearAlgebra;

for dimension = 2:3
    # setup
    mesh = []
    if dimension == 2
        mesh = Mesh2D(1.0, 1.0)
    elseif dimension == 3
        mesh = Mesh3D(1.0, 1.0, 1.0)
    end
    diffusion = GalleryOperator("diffusion", 3, 3, mesh)
    bddc = DirichletBDDC(diffusion)

    # compute symbols
    A = computesymbols(bddc, [0.2], π * ones(dimension))

    # verify
    eigenvalues = real(eigvals(A))
    if dimension == 2
        @assert minimum(eigenvalues) ≈ 0.7999999999999998
        @assert maximum(eigenvalues) ≈ 0.8
    elseif dimension == 3
        @assert minimum(eigenvalues) ≈ 0.7801226993865031
        @assert maximum(eigenvalues) ≈ 0.8
    end
end

# output
source

Dirichlet BDDC

Examples

This is an example of a simple Dirichlet BDDC preconditioner.

# ------------------------------------------------------------------------------
# Dirichlet BDDC example
# ------------------------------------------------------------------------------

using LFAToolkit
using LinearAlgebra

# setup
mesh = Mesh2D(1.0, 1.0)
p = 3

# diffusion operator
diffusion = GalleryOperator("diffusion", p + 1, p + 1, mesh)

# Dirichlet BDDC preconditioner
bddc = DirichletBDDC(diffusion)

# compute operator symbols
A = computesymbols(bddc, [0.2], [π, π])
eigenvalues = real(eigvals(A))

# ------------------------------------------------------------------------------

This is an example of a simple Dirichlet BDDC preconditioner on a macro-element patch.

# ------------------------------------------------------------------------------
# lumped BDDC on macro elements example
# ------------------------------------------------------------------------------

using LFAToolkit
using LinearAlgebra

# setup
mesh = Mesh2D(1.0, 1.0)
p = 1
numberelements1d = 4

# operator
diffusion = GalleryMacroElementOperator("diffusion", p + 1, p + 2, numberelements1d, mesh)

# BDDC smoother
bddc = DirichletBDDC(diffusion)

# compute operator symbols
A = computesymbols(bddc, [0.2], [π, π])
eigenvalues = real(eigvals(A))

# ------------------------------------------------------------------------------

Documentation

LFAToolkit.DirichletBDDCFunction
DirichletBDDC(operator)

Dirichlet BDDC preconditioner for finite element operators

Arguments:

  • operator::Operator: finite element operator to precondition

Returns:

  • Dirichlet BDDC preconditioner object

Example:

# setup
mesh = Mesh2D(1.0, 1.0);

# operator
finediffusion = GalleryOperator("diffusion", 5, 5, mesh);

# preconditioner
bddc = DirichletBDDC(finediffusion);

# verify
println(bddc)
println(bddc.operator)

# output

Dirichlet BDDC preconditioner
finite element operator:
2d mesh:
    dx: 1.0
    dy: 1.0

2 inputs:
operator field:
  tensor product basis:
    numbernodes1d: 5
    numberquadraturepoints1d: 5
    numbercomponents: 1
    dimension: 2
  evaluation mode:
    gradient
operator field:
  tensor product basis:
    numbernodes1d: 5
    numberquadraturepoints1d: 5
    numbercomponents: 1
    dimension: 2
  evaluation mode:
    quadratureweights

1 output:
operator field:
  tensor product basis:
    numbernodes1d: 5
    numberquadraturepoints1d: 5
    numbercomponents: 1
    dimension: 2
  evaluation mode:
    gradient
source