Skip to content

Assemblers

All docstrings from Mantis.Assemblers

Mantis.Assemblers Module
julia
module Assemblers

Contains all assembly-related structs and functions.

source
Mantis.Assemblers.WeakForm Type
julia
WeakForm{manifold_dim, LHS, RHS, I}

Structure representing the weak-formulation of a continuous Petrov-Galerking method. Both the left and right-hand sides of the formulation are given as blocks of real-valued operators; these are defined from a set of inputs holding the test, trial and forcing terms, and a constructor method that defines where the blocks are placed.

Fields

  • lhs_expressions::LHS: The left-hand side blocks of the weak-formulation.

  • rhs_expressions::RHS: The right-hand side blocks of the weak-formulation.

  • inputs::I: The inputs for the weak-formulation, which include the test and trial spaces, and forcing terms.

Type parameters

  • manifold_dim::Int: The dimension of the manifold where the weak-formulation is defined.

  • LHS: The type of the left-hand side expressions. Each row-column entry should be a subtype of AbstractRealValuedOperator or 0.

  • RHS: The type of the right-hand side expressions. Each row-column entry should be a subtype of AbstractRealValuedOperator or 0.

  • I: The type of the inputs. It should be a subtype of WeakFormInputs{manifold_dim}.

Inner constructors

  • WeakForm(inputs::I, constructor::F): Creates a new WeakForm instance with the given inputs and constructor function. The constructor function is used to generate the left-hand side and right-hand side blocks of real-valued operators.
source
Mantis.Assemblers.WeakFormInputs Type
julia
WeakFormInputs{manifold_dim, TeF, TrF, F} <: AbstractInputs

Container for test and trial spaces, and forcing terms to be used in a weak-formulation.

Fields

  • test_forms::TeF: The test forms for the weak-formulation.

  • trial_forms::TrF: The trial forms for the weak-formulation.

  • forcings::F: The forcing terms for the weak-formulation, possibly nothing.

Type parameters

  • manifold_dim::Int: The dimension of the manifold where the weak-formulation is defined.

  • TeF: The type of the tupe of test forms. Each entry should be a subtype of Forms.AbstractFormSpace.

  • TrF: The type of the tupe of trial forms. Each entry should be a subtype of Forms.AbstractFormSpace.

  • F: The type of the tupe of forcing terms. Each entry should be a subtype of Forms.AbstractFormField.

Inner constructors

  • WeakFormInputs(test_forms::TeF, trial_forms::TrF, forcings::F): Creates a new WeakFormInputs instance with the given test forms, trial forms, and forcing terms.

  • WeakFormInputs(test_forms::TeF, trial_forms::TrF): Creates a new WeakFormInputs instance with the given test forms and trial forms. The forcing terms are set to nothing.

  • WeakFormInputs(test_forms::TeF, trial_forms::TrF, forcing::F): Creates a new WeakFormInputs instance from a single test and trial space and a forcing term.

  • WeakFormInputs(test_forms::TeF, trial_forms::TrF): Creates a new WeakFormInputs instance from a single test and trial space. The forcing terms are set to nothing.

  • WeakFormInputs(forms::TrF, forcing::F): Creates a new WeakFormInputs instance with the given trial spaces and forcing terms. The test spaces are set to the same as the trial spaces.

  • WeakFormInputs(forms::TrF): Creates a new WeakFormInputs instance with the given trial spaces. The test spaces are set to the same as the trial spaces and the forcing terms are set to nothing.

  • WeakFormInputs(forms::TrF, forcing::F): Creates a new WeakFormInputs instance with a single trial space and forcing term. The test space is set to the same as the trial space.

  • WeakFormInputs(forms::TrF): Creates a new WeakFormInputs instance with a single trial space. The test space is set to the same as the trial space and the forcing term is set to nothing.

source
Mantis.Assemblers.L2_projection Method
julia
L2_projection(inputs::AbstractInputs, dΩ::Quadrature.AbstractGlobalQuadratureRule)

Function to compute the L2 projection of a function onto a discrete form space.

Arguments

  • inputs::AbstractInputs: The inputs for the weak form assembly, including test, trial and forcing terms.

  • dΩ::Quadrature.AbstractGlobalQuadratureRule: The quadrature rule to use for the integral evaluation.

Returns

  • lhs_expression<:NTuple{num_lhs_rows, NTuple{num_lhs_cols, AbstractRealValuedOperator}}: The left-hand side of the weak form, which is a tuple of tuples contain all the blocks of the left-hand side matrix.

  • rhs_expression<:NTuple{num_rhs_rows, NTuple{num_rhs_cols, AbstractRealValuedOperator}}: The right-hand side of the weak form, which is a tuple of tuples contain all the blocks of the right-hand side matrix.

source
Mantis.Assemblers.add_bc! Method
julia
add_bc!(
    lhs_rows::Vector{Int},
    lhs_cols::Vector{Int},
    lhs_vals::Vector{Float64},
    rhs_rows::Vector{Int},
    rhs_cols::Vector{Int},
    rhs_vals::Vector{Float64},
    dirichlet_bcs::Dict{Int, Float64},
)

Adds Dirichlet boundary conditions to the left-hand side and right-hand side matrices.

Arguments

  • lhs_rows::Vector{Int}: The row indices of the left-hand side matrix.

  • lhs_cols::Vector{Int}: The column indices of the left-hand side matrix.

  • lhs_vals::Vector{Float64}: The values of the left-hand side matrix.

  • rhs_rows::Vector{Int}: The row indices of the right-hand side matrix.

  • rhs_cols::Vector{Int}: The column indices of the right-hand side matrix.

  • rhs_vals::Vector{Float64}: The values of the right-hand side matrix.

  • dirichlet_bcs::Dict{Int, Float64}: The Dirichlet boundary conditions, where the key is the index of the boundary condition and the value is the boundary condition value.

Returns

  • lhs_rows::Vector{Int}: The updated row indices of the left-hand side matrix.

  • lhs_cols::Vector{Int}: The updated column indices of the left-hand side matrix.

  • lhs_vals::Vector{Float64}: The updated values of the left-hand side matrix.

  • rhs_rows::Vector{Int}: The updated row indices of the right-hand side matrix.

  • rhs_cols::Vector{Int}: The updated column indices of the right-hand side matrix.

  • rhs_vals::Vector{Float64}: The updated values of the right-hand side matrix.

source
Mantis.Assemblers.add_expression_contributions! Method
julia
add_expression_contributions!(
    rows::Vector{Int},
    cols::Vector{Int},
    vals::Vector{Float64},
    counts::Int,
    expression,
    element_id::Int,
    test_offset::Int,
    trial_offset::Int,
)

Updates the row, column, and value vectors with contributions from the specified real-valued expression at the element given by element_id.

Arguments

  • rows::Vector{Int}: The row indices of the matrix.

  • cols::Vector{Int}: The column indices of the matrix.

  • vals::Vector{Float64}: The values of the matrix.

  • counts::Int: The current count of non-zero entries in the matrix.

  • expressions: The expression to evaluate.

  • element_id::Int: The identifier of the element.

  • test_offsets::Int: The offset for the test function.

  • trial_offsets::Int: The offset for the trial function.

Returns

  • rows::Vector{Int}: The updated row indices of the matrix.

  • cols::Vector{Int}: The updated column indices of the matrix.

  • vals::Vector{Float64}: The updated values of the matrix.

  • counts::Int: The updated count of non-zero entries in the matrix.

source
Mantis.Assemblers.analytical_maxwell_eigenfunction Method
julia
analytical_maxwell_eigenfunction(
    m::Int, n::Int, scale_factors::NTuple{2, Float64}, x::Matrix{Float64}
)

Evaluates the analytical Maxwell eigenfunction for the eigenmode (m, n) at points x.

Arguments

  • m::Int: Eigenmode component x.

  • n::Int: Eigenmode component y.

  • scale_factors::NTuple{2, Int}: Scaling factors based on the size of the domain.

  • x::Matrix{Float64}: Evaluation points.

Returns

  • NTuple{2, Vector{Float64}}: The evaluated Maxwell eigenfunction at points x.
source
Mantis.Assemblers.assemble Method
julia
assemble(
    weak_form::WeakForm{manifold_dim, LHS, RHS, I},
    dirichlet_bcs::Dict{Int, Float64}=Dict{Int, Float64}();
    lhs_type::Type=spa.SparseMatrixCSC{Float64, Int},
    rhs_type::Type=Matrix{Float64},
) where {
    manifold_dim,
    num_rows,
    lhs_num_cols,
    rhs_num_cols,
    LHS <: NTuple{
        num_rows, NTuple{lhs_num_cols, Union{Int, Forms.AbstractRealValuedOperator}}
    },
    RHS <: NTuple{
        num_rows, NTuple{rhs_num_cols, Union{Int, Forms.AbstractRealValuedOperator}}
    },
    I,
}

Assemble the left- and right-hand sides of a discrete Petrov-Galerkin problem for the given weak-formulation and Dirichlet boundary conditions.

Arguments

  • weak_form::WeakForm{manifold_dim, LHS, RHS, I}: The weak form to assemble.

  • dirichlet_bcs::Dict{Int, Float64}: A dictionary containing the Dirichlet boundary conditions, where the key is the index of the boundary condition and the value is the boundary condition value.

  • lhs_type::Type: The type of the left-hand side matrix. Default is SparseMatrixCSC{Float64, Int}.

  • rhs_type::Type: The type of the right-hand side matrix. Default is Matrix{Float64}.

Returns

  • lhs::lhs_type: The assembled left-hand side matrix.

  • rhs::rhs_type: The assembled right-hand side vector.

source
Mantis.Assemblers.build_matrix Method
julia
build_matrix(
    matrix_type::Type,
    rows::Vector{Int},
    cols::Vector{Int},
    vals::Vector{Float64},
    size::Tuple{Int, Int},
)

Returns a matrix of the specified type with the given row and column indices and values.

Arguments

  • matrix_type::Type: The type of matrix to build.

  • rows::Vector{Int}: The row indices of the matrix.

  • cols::Vector{Int}: The column indices of the matrix.

  • vals::Vector{Float64}: The values of the matrix.

  • size::Tuple{Int, Int}: The size of the matrix.

Returns

  • ::matrix_type: The constructed matrix of the specified type.
source
Mantis.Assemblers.get_analytical_maxwell_eig Method
julia
get_maxwell_eig(num_eig::Int, geom::G) where {G <: Geometry.AbstractGeometry{2}}

Returns the first num_eig eigenvalues and 1-form eigenfunctions on the geometry geom.

Arguments

  • num_eig::Int: Number of eigenvalues and eigenfunctions to compute.

  • geom::Geometry.AbstractGeometry{2}: The two-dimensional geometry.

Returns

  • Vector{Float64}: The first num_eig analytical eigenvalues.

  • Vector{Forms.AnalyticalFormField{2, 1, G}: The first num_eig analytical eigenfunctions.

source
Mantis.Assemblers.get_estimated_nnz_per_elem Method
julia
get_estimated_nnz_per_elem(wf::WeakForm)

Returns the estimated number of non-zero entries per element for the left- and right-hand sides of the weak-formulation.

Arguments

  • wf::WeakForm: The weak-formulation for which the estimated number of non-zero entries is to be determined.

Returns

  • Tuple(Int, Int): The estimated number of non-zero entries per element for the left-hand side and right-hand side of the weak-formulation, respectively.
source
Mantis.Assemblers.get_num_elements Method
julia
get_num_elements(wf::WeakForm)

Returns the number of elements over which the discrete weak-formulation is defined.

Arguments

  • wf::WeakForm: The weak-formulation for which the number of elements is to be determined.

Returns

  • ::Int: The number of elements.
source
Mantis.Assemblers.get_num_evaluation_elements Method
julia
get_num_evaluation_elements(wf::WeakForm)

Returns the maximum number of elements over which the weak form blocks are evaluated. This is the max over all lhs and rhs expression blocks.

Arguments

  • wf::WeakForm: The weak-formulation for which the number of quadrature elements is to be determined.

Returns

  • ::Int: The number of quadrature elements.
source
Mantis.Assemblers.get_pre_allocation Method
julia
get_pre_allocation(weak_form::WeakForm, side::String)

Returns pre-allocated row, column, and value vectors for the left-hand side (lhs) or right-hand side (rhs) matrix.

Arguments

  • weak_form::WeakForm: The weak form to use for the pre-allocation.

  • side::String: The side of the matrix to pre-allocate. Must be either "lhs" or "rhs".

Returns

  • rows::Vector{Int}: The pre-allocated row indices.

  • cols::Vector{Int}: The pre-allocated column indices.

  • vals::Vector{Float64}: The pre-allocated values.

source
Mantis.Assemblers.get_test_offsets Method
julia
get_test_offsets(
    wf::WeakForm{manifold_dim, LHS, RHS, I}
) where {
    manifold_dim,
    num_rows,
    lhs_num_cols,
    rhs_num_cols,
    LHS <:
    NTuple{
        num_rows, NTuple{lhs_num_cols, Union{Int, Forms.AbstractRealValuedOperator}}
    },
    RHS <:
    NTuple{
        num_rows, NTuple{rhs_num_cols, Union{Int, Forms.AbstractRealValuedOperator}}
    },
    I,
}

Returns the offsets of the test function spaces used in the weak form.

Arguments

  • wf::WeakForm{manifold_dim, LHS, RHS, I}: The weak form being used.

Returns

  • NTuple{num_rows, Int}: The offsets of the test function spaces.
source
Mantis.Assemblers.get_trial_offsets Method
julia
get_trial_offsets(
    wf::WeakForm{manifold_dim, LHS, RHS, I}
) where {
    manifold_dim,
    num_rows,
    lhs_num_cols,
    rhs_num_cols,
    LHS <:
    NTuple{
        num_rows, NTuple{lhs_num_cols, Union{Int, Forms.AbstractRealValuedOperator}}
    },
    RHS <:
    NTuple{
        num_rows, NTuple{rhs_num_cols, Union{Int, Forms.AbstractRealValuedOperator}}
    },
    I,
}

Returns the offsets of the trial function spaces used in the weak form.

Arguments

  • wf::WeakForm{manifold_dim, LHS, RHS, I}: The weak form being used.

Returns

  • NTuple{lhs_num_cols, lhs_num_cols}: The offsets of the trial function spaces.
source
Mantis.Assemblers.maxwell_eigenvalue Method
julia
maxwell_eigenvalue(inputs::WeakFormInputs, dΩ::Quadrature.AbstractGlobalQuadratureRule)

Function for assembling the weak form of the Maxwell eigenvalue problem.

Arguments

  • inputs::WeakFormInputs: The inputs for the weak form assembly, including test and trial spaces.

  • dΩ::Quadrature.AbstractGlobalQuadratureRule: The quadrature rule to use for the integral evaluation.

Returns

  • lhs_expression<:NTuple{num_lhs_rows, NTuple{num_lhs_cols, AbstractRealValuedOperator}}: The left-hand side of the weak form, which is a tuple of tuples contain all the blocks of the left-hand side matrix.

  • rhs_expression<:NTuple{num_rhs_rows, NTuple{num_rhs_cols, AbstractRealValuedOperator}}: The right-hand side of the weak form, which is a tuple of tuples contain all the blocks of the right-hand side matrix.

source
Mantis.Assemblers.n_form_hodge_laplacian Method
julia
n_form_hodge_laplacian(
    inputs::WeakFormInputs, dΩ::Quadrature.AbstractGlobalQuadratureRule
)

Function for assembling the weak form of the n-form Hodge Laplacian problem.

Arguments

  • inputs::WeakFormInputs: The inputs for the weak form assembly, including test, trial and forcing terms.

  • dΩ::Quadrature.AbstractGlobalQuadratureRule: The quadrature rule to use for the integral evaluation.

Returns

  • lhs_expressions<:NTuple{num_lhs_rows, NTuple{num_lhs_cols, AbstractRealValuedOperator}}: The left-hand side of the weak form, which is a tuple of tuples contain all the blocks of the left-hand side matrix.

  • rhs_expressions<:NTuple{num_rhs_rows, NTuple{num_rhs_cols, AbstractRealValuedOperator}}: The right-hand side of the weak form, which is a tuple of tuples contain all the blocks of the right-hand side matrix.

source
Mantis.Assemblers.one_form_hodge_laplacian Method
julia
one_form_hodge_laplacian(
    inputs::AbstractInputs, dΩ::Quadrature.AbstractGlobalQuadratureRule
)

Function for assembling the weak form of the 1-form Hodge Laplacian problem.

Arguments

  • inputs::AbstractInputs: The inputs for the weak form assembly, including test, trial and forcing terms.

  • dΩ::Quadrature.AbstractGlobalQuadratureRule: The quadrature rule to use for the integral evaluation.

Returns

  • lhs_expressions<:NTuple{num_lhs_rows, NTuple{num_lhs_cols, AbstractRealValuedOperator}}: The left-hand side of the weak form, which is a tuple of tuples contain all the blocks of the left-hand side matrix.

  • rhs_expressions<:NTuple{num_rhs_rows, NTuple{num_rhs_cols, AbstractRealValuedOperator}}: The right-hand side of the weak form, which is a tuple of tuples contain all the blocks of the right-hand side matrix.

source
Mantis.Assemblers.solve_L2_projection Method
julia
solve_L2_projection(Xᵏ, fₑ, dΩ)

Returns the solution of the weak form of the L2 projection.

Arguments

  • Xᵏ: The k-form space to use as trial and test space.

  • fₑ: The forcing term to use for the right-hand side of the weak formulation.

  • : The quadrature rule to use for the assembly.

Returns

  • fₕ::FormField: The projection of fₑ onto Xᵏ.
source
Mantis.Assemblers.solve_maxwell_eig Method
julia
solve_maxwell_eig(
    X⁰::Forms.AbstractFormSpace{2, 0},
::Forms.AbstractFormSpace{2, 1},
::Quadrature.AbstractGlobalQuadratureRule{2},
    num_eig::Int;
    verbose::Bool=false,
)

Returns the first num_eig eigenvalues and 1-form eigenfunctions of the Maxwell eigenvalue problem.

Arguments

  • X⁰::Forms.AbstractFormSpace{2, 0}: The 0-form space to use as trial and test space.

  • X¹::Forms.AbstractFormSpace{2, 1}: The 1-form space to use as trial and test space.

  • dΩ::Quadrature.AbstractGlobalQuadratureRule{2}: The quadrature rule to use for the assembly.

  • num_eig::Int: The number of eigenvalues and eigenfunctions to compute.

  • verbose::Bool=false: Whether to print the nullspace offset.

Returns

  • ω²ₕ::Vector{Float64}: The first num_eig eigenvalues.

  • u¹ₕ::Vector{Forms.FormField{2, 1}}: The first num_eig eigenfunctions.

source
Mantis.Assemblers.solve_maxwell_eig Method
julia
solve_maxwell_eig(
    complex::C,
    dΩₐ::Quadrature.StandardQuadrature{manifold_dim},
    num_steps::Int,
    dorfler_parameter::Float64,
    dΩₑ::Quadrature.StandardQuadrature{manifold_dim},
    Lchains::Bool,
    eigenfunction::Int,
    num_eig::Int,
    scale_factors::NTuple{manifold_dim, Float64};
    verbose::Bool=false,
) where {manifold_dim, num_forms, C <: NTuple{num_forms, Forms.AbstractFormSpace}}

Returns the first num_eig eigenvalues and 1-form eigenfunctions of the Maxwell eigenvalue problem for an adaptive loop.

Arguments

  • complex::C: The initial de Rham complex to use for the problem.

  • dΩₐ::Quadrature.StandardQuadrature{manifold_dim}: The quadrature rule to use for the assembly.

  • num_steps::Int: The number of adaptive steps to perform.

  • dorfler_parameter::Float64: The Dörfler marking parameter.

  • dΩₑ::Quadrature.StandardQuadrature{manifold_dim}: The quadrature rule to use for the error estimation.

  • Lchains::Bool: Whether to use L-chains for the refinement.

  • eigenfunction::Int: The index of the eigenfunction to use for the error estimation.

  • num_eig::Int: The number of eigenvalues and eigenfunctions to compute.

  • scale_factors::NTuple{manifold_dim, Float64}: The scaling factors for the geometry.

  • verbose::Bool=false: Whether to print the progress of the adaptive loop.

Returns

  • ω²ₕ::Vector{Float64}: The first num_eig eigenvalues.

  • u¹ₕ::Vector{Forms.FormField{manifold_dim, 1, G}}: The first num_eig eigenfunctions.

source
Mantis.Assemblers.solve_one_form_hodge_laplacian Function
julia
solve_one_form_hodge_laplacian(X⁰, X¹, f¹, dΩ)

Returns the solution of the weak form of the 1-form Hodge Laplacian.

Arguments

  • X⁰: The 0-form space to use as trial and test space.

  • : The 1-form space to use as trial and test space.

  • : The forcing term to use for the right-hand side of the weak formulation.

  • : The quadrature rule to use for the assembly.

Returns

  • δu¹ₕ::Forms.FormField: The 0-form solution of the weak-formulation.

  • u¹ₕ::Forms.FormField: The 1-form solution of the weak-formulation.

source
Mantis.Assemblers.solve_one_form_hodge_laplacian Method
julia
solve_one_form_hodge_laplacian(
    complex::C,
    forcing_function::Function,
    dΩₐ::Quadrature.AbstractGlobalQuadratureRule{manifold_dim},
    num_steps::Int,
    dorfler_parameter::Float64,
    dΩₑ::Quadrature.AbstractGlobalQuadratureRule{manifold_dim},
    Lchains::Bool;
    verbose::Bool=false
) where {manifold_dim, num_forms, C <: NTuple{num_forms, Forms.AbstractFormSpace}}

Returns the solution of the weak form of the 1-form Hodge Laplacian from an adaptive loop.

Arguments

  • complex::C: The initial de Rham complex to use for the problem.

  • forcing_function::Function: The function to use for the forcing term.

  • dΩₐ::Quadrature.AbstractGlobalQuadratureRule{manifold_dim}: The quadrature rule to use for the assembly.

  • num_steps::Int: The number of steps to use for the adaptive loop.

  • dorfler_parameter::Float64: The parameter to use for the Dörfler marking.

  • dΩₑ::Quadrature.AbstractGlobalQuadratureRule{manifold_dim}: The quadrature rule to use for the error estimation.

  • Lchains::Bool: Whether to use L-chains for the refinement.

  • verbose::Bool=false: Whether to print the progress of the adaptive loop.

source
Mantis.Assemblers.solve_volume_form_hodge_laplacian Method
julia
solve_volume_form_hodge_laplacian(Xⁿ⁻¹, Xⁿ, fₑ, dΩ)

Returns the solution of the weak form of the n-form Hodge Laplacian.

Arguments

  • Xⁿ⁻¹: The (n-1)-form space to use as trial and test space.

  • Xⁿ: The n-form space to use as trial and test space.

  • fₑ: The forcing term to use for the right-hand side of the weak formulation.

  • : The quadrature rule to use for the assembly.

Returns

  • u¹ₕ: The (n-1)-form solution of the weak-formulation.

  • ϕ²ₕ: The n-form solution of the weak-formulation.

source
Mantis.Assemblers.solve_zero_form_hodge_laplacian Method
julia
solve_zero_form_hodge_laplacian(X⁰, fₑ, dΩ)

Returns the solution of the weak form of the 0-form Hodge Laplacian.

Arguments

  • X⁰: The 0-form space to use as trial and test space.

  • fₑ: The forcing term to use for the right-hand side of the weak formulation.

  • : The quadrature rule to use for the assembly.

Returns

  • ::Forms.FormField: The solution of the weak-formulation.
source
Mantis.Assemblers.zero_form_hodge_laplacian Method
julia
zero_form_hodge_laplacian(
    inputs::AbstractInputs, dΩ::Quadrature.AbstractGlobalQuadratureRule
)

Function for assembling the weak form of the 0-form Hodge Laplacian.

Arguments

  • inputs::AbstractInputs: The inputs for the weak form assembly, including test, trial and forcing terms.

  • dΩ::Quadrature.AbstractGlobalQuadratureRule: The quadrature rule to use for the integral evaluation.

Returns

  • lhs_expression<:NTuple{num_lhs_rows, NTuple{num_lhs_cols, AbstractRealValuedOperator}}: The left-hand side of the weak form, which is a tuple of tuples contain all the blocks of the left-hand side matrix.

  • rhs_expression<:NTuple{num_rhs_rows, NTuple{num_rhs_cols, AbstractRealValuedOperator}}: The right-hand side of the weak form, which is a tuple of tuples contain all the blocks of the right-hand side matrix.

source