Assemblers
All docstrings from Mantis.Assemblers
Mantis.Assemblers Module
module AssemblersContains all assembly-related structs and functions.
sourceMantis.Assemblers.WeakForm Type
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 ofAbstractRealValuedOperatoror0.RHS: The type of the right-hand side expressions. Each row-column entry should be a subtype ofAbstractRealValuedOperatoror0.I: The type of the inputs. It should be a subtype ofWeakFormInputs{manifold_dim}.
Inner constructors
WeakForm(inputs::I, constructor::F): Creates a newWeakForminstance 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.
Mantis.Assemblers.WeakFormInputs Type
WeakFormInputs{manifold_dim, TeF, TrF, F} <: AbstractInputsContainer 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 ofForms.AbstractFormSpace.TrF: The type of the tupe of trial forms. Each entry should be a subtype ofForms.AbstractFormSpace.F: The type of the tupe of forcing terms. Each entry should be a subtype ofForms.AbstractFormField.
Inner constructors
WeakFormInputs(test_forms::TeF, trial_forms::TrF, forcings::F): Creates a newWeakFormInputsinstance with the given test forms, trial forms, and forcing terms.WeakFormInputs(test_forms::TeF, trial_forms::TrF): Creates a newWeakFormInputsinstance 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 newWeakFormInputsinstance from a single test and trial space and a forcing term.WeakFormInputs(test_forms::TeF, trial_forms::TrF): Creates a newWeakFormInputsinstance from a single test and trial space. The forcing terms are set to nothing.WeakFormInputs(forms::TrF, forcing::F): Creates a newWeakFormInputsinstance 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 newWeakFormInputsinstance 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 newWeakFormInputsinstance 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 newWeakFormInputsinstance 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.
Mantis.Assemblers.L2_projection Method
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.
Mantis.Assemblers.add_bc! Method
add_bc!(lhs::AbstractMatrix, rhs::AbstractMatrix, dirichlet_bcs::Dict)Remove the rows and columns of lhs and rhs as specified by the keys of dirichlet_bcs. See also set_diagonal!. ```
Mantis.Assemblers.add_bc! Method
add_bc!(
lhs::AbstractMatrix{T}, rhs::AbstractVector{T}, dirichlet_bcs::Dict{Int, T}
) where {T}Set the rows to 0 and diagonals to 1 in lhs, and rows of rhs to a value, as specified by the row indices (keys) and values of dirichlet_bcs.
Examples
using Mantis
Assemblers.add_bc!([2.0 2.0 2.0; 2.0 2. 2.0; 2.0 2.0 2.0], zeros(3), Dict(2 => 42.0))
# output
([2.0 2.0 2.0; 0.0 1.0 0.0; 2.0 2.0 2.0], [0.0, 42.0, 0.0])Mantis.Assemblers.add_expression_contributions! Method
add_expression_contributions!(
rows::Vector{Int},
cols::Vector{Int},
vals::AbstractVector,
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 array.cols::Vector{Int}: The column indices of the array.vals::AbstractVector: The values of the array.counts::Int: The current count of non-zero entries in the array.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 array.cols::Vector{Int}: The updated column indices of the array.vals::AbstractVector: The updated values of the array.counts::Int: The updated count of non-zero entries in the array.
Mantis.Assemblers.analytical_maxwell_eigenfunction Method
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 pointsx.
Mantis.Assemblers.assemble Method
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=Vector{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 array. Default isSparseMatrixCSC{Float64, Int}.rhs_type::Type: The type of the right-hand side array. Default isVector{Float64}.
Returns
lhs::lhs_type: The assembled left-hand side array.rhs::rhs_type: The assembled right-hand side vector.
Mantis.Assemblers.build_array Method
build_array(
array_type::Type{A},
rows::Vector{Int},
cols::Vector{Int},
vals::AbstractVector,
size::Tuple{Int, Int},
) where {A <: AbstractArray}Returns a array of the specified type with the given row and column indices and values.
Arguments
array_type::Type{AbstractArray}: The type of array to build.rows::Vector{Int}: The row indices of the array.cols::Vector{Int}: The column indices of the array.vals::AbstractVector: The values of the array.size::Tuple{Int, Int}: The size of the array.
Returns
::array_type: The constructed array of the specified type.
Mantis.Assemblers.get_analytical_maxwell_eig Method
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 firstnum_eiganalytical eigenvalues.Vector{Forms.AnalyticalFormField{2, 1, G}: The firstnum_eiganalytical eigenfunctions.
Mantis.Assemblers.get_estimated_nnz_per_elem Method
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.
Mantis.Assemblers.get_num_elements Method
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.
Mantis.Assemblers.get_num_evaluation_elements Method
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.
Mantis.Assemblers.get_pre_allocation Method
get_pre_allocation(
weak_form::WeakForm, side::String, ::Type{A}
) where {T, A <: AbstractArray{T}}Returns pre-allocated row, column, and value vectors for the left-hand side (lhs) or right-hand side (rhs) array.
Arguments
weak_form::WeakForm: The weak form to use for the pre-allocation.side::String: The side of the array 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{T}: The pre-allocated values.
Mantis.Assemblers.get_test_offsets Method
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.
Mantis.Assemblers.get_trial_offsets Method
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.
Mantis.Assemblers.maxwell_eigenvalue Method
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.
Mantis.Assemblers.n_form_hodge_laplacian Method
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.
Mantis.Assemblers.one_form_hodge_laplacian Method
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.
Mantis.Assemblers.set_diagonal! Method
set_diagonal!(lhs::AbstractMatrix, rhs::AbstractMatrix, dirichlet_bcs::Dict)Remove the rows and columns of lhs and rhs as specified by the keys of dirichlet_bcs.
Examples
using Mantis
Assemblers.set_diagonal!(
[2. 0. 0.; 0. 2. 0.; 0. 0. 2.], [2. 0. 0.; 0. 2. 0.; 0. 0. 2.], Dict(2 => 42.0)
)
# output
([2.0 0.0; 0.0 2.0], [2.0 0.0; 0.0 2.0])Mantis.Assemblers.set_diagonal! Method
set_diagonal!(
lhs::AbstractMatrix{T}, rhs::AbstractVector{T}, dirichlet_bcs::Dict{Int, T}
) where {T}Set the diagonals of lhs to 1 and rows of rhs to a value, as specified by the row indices (keys) and values of dirichlet_bcs. See also add_bc!.
Examples
using Mantis
Assemblers.set_diagonal!([2.0 2.0 2.0; 2.0 2. 2.0; 2.0 2.0 2.0], zeros(3), Dict(2 => 42.0))
# output
([2.0 2.0 2.0; 2.0 1.0 2.0; 2.0 2.0 2.0], [0.0, 42.0, 0.0])Mantis.Assemblers.solve_L2_projection Method
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.dΩ: The quadrature rule to use for the assembly.
Returns
fₕ::FormField: The projection offₑontoXᵏ.
Mantis.Assemblers.solve_maxwell_eig Method
solve_maxwell_eig(
X⁰::Forms.AbstractFormSpace{2, 0},
X¹::Forms.AbstractFormSpace{2, 1},
dΩ::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 firstnum_eigeigenvalues.u¹ₕ::Vector{Forms.FormField{2, 1}}: The firstnum_eigeigenfunctions.
Mantis.Assemblers.solve_maxwell_eig Method
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 firstnum_eigeigenvalues.u¹ₕ::Vector{Forms.FormField{manifold_dim, 1, G}}: The firstnum_eigeigenfunctions.
Mantis.Assemblers.solve_one_form_hodge_laplacian Function
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.X¹: The 1-form space to use as trial and test space.f¹: The forcing term to use for the right-hand side of the weak formulation.dΩ: 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.
Mantis.Assemblers.solve_one_form_hodge_laplacian Method
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.
Mantis.Assemblers.solve_volume_form_hodge_laplacian Method
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.dΩ: 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.
Mantis.Assemblers.solve_zero_form_hodge_laplacian Method
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.dΩ: The quadrature rule to use for the assembly.
Returns
::Forms.FormField: The solution of the weak-formulation.
Mantis.Assemblers.zero_form_hodge_laplacian Method
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.
Mantis.Assemblers.zero_rows! Method
zero_rows!(
lhs_vals::AbstractVector{T},
rhs_vals::AbstractVector{T},
lhs_rows::Vector{Int},
rhs_rows::Vector{Int},
dirichlet_bcs::Dict{Int, T},
) where {T}For each index i that is a key of dirichlet_bcs, set the corresponding values from that row to zero, both in lhs_vals and rhs_vals.
Examples
using Mantis
Assemblers.zero_rows!([1., 1., 1.], [1., 2., 3.], [1, 2, 3], [1, 3, 2], Dict(2 => 42.0))
# output
([1.0, 0.0, 1.0], [1.0, 2.0, 0.0])