Assembly namespace
Module for performing the assembly of linear algebra objects from variational expressions.
The Assembly module provides the core infrastructure for assembling finite element matrices and vectors from variational forms. It supports both sequential and parallel assembly strategies, with specializations for different operator types and finite element spaces.
Assembly Strategies
| Strategy | Class | Description |
|---|---|---|
| Sequential | Assembly:: | Single-threaded assembly. Iterates over mesh cells, evaluates local element matrices/vectors, and accumulates into the global system. Deterministic and reproducible. |
| OpenMP | Assembly::OpenMP | Shared-memory parallel assembly using OpenMP. Each thread processes a subset of cells and accumulates contributions using thread-safe operations. Requires compilation with RODIN_USE_OPENMP=ON. |
| Default | Assembly:: | Compile-time policy selector. Maps to OpenMP when OpenMP is available, otherwise falls back to Sequential. Users rarely need to specify this directly. |
What Gets Assembled
- Bilinear Form : Assembles into a sparse matrix where are trial basis functions and are test basis functions.
- Linear Form : Assembles into a vector .
- Dirichlet BC: Assembles into an
IndexMapmapping constrained DOFs to their prescribed values. - Problem: Orchestrates the full assembly pipeline — bilinear form, linear form(s), Dirichlet elimination — into a complete LinearSystem.
Assembly Pipeline
Assembly is normally triggered automatically when Problem::solve() is called. The pipeline is:
- Assemble the bilinear form into a sparse matrix
- Assemble the linear form into a vector
- Apply Dirichlet boundary conditions (row/column elimination)
- Package into a
LinearSystemand pass to the solver
Users can also trigger assembly manually by calling BilinearForm::assemble() or LinearForm::assemble().
Multi-Field Assembly
For saddle-point problems (e.g. Stokes), the Assembly module supports block-structured systems via Tuple<BilinearForm...> and Tuple<LinearForm...>, which assemble into block matrices and vectors. This is used by multi-field Problem objects with multiple trial/test function pairs.
Classes
-
template<class LinearAlgebraType, class Operand>class AssemblyBase
- Base class template for assembly operations.
-
template<class OperatorType, class Solution, class TrialFES, class TestFES>class AssemblyBase<OperatorType, Variational::BilinearForm<Solution, TrialFES, TestFES, OperatorType>>
- Base class for bilinear form assembly operations.
-
template<class VectorType, class FES>class AssemblyBase<VectorType, Variational::LinearForm<FES, VectorType>>
- Base class for linear form assembly operations.
-
template<class TrialFES, class TestFES>class BilinearFormAssemblyInput
- Input data for bilinear form assembly.
-
template<class ... Ts>class BilinearFormTupleAssemblyInput
- Input data for tuple of bilinear forms assembly.
-
template<class ... Ts>class Default
- Default assembly strategy selector.
-
template<>class Default<Context::Local>
- Default assembly strategy for local context without OpenMP.
-
template<>class Default<Context::Local, Context::Local>
- Default assembly strategy for mixed local contexts without OpenMP.
-
template<>class Default<Context::MPI>
- Selects MPI assembly when only the trial/test context is MPI.
-
template<>class Default<Context::MPI, Context::MPI>
- Selects MPI assembly when both operand and assembly contexts are MPI.
-
template<class Scalar, class Solution, class FES, class Value>class DirichletBCAssemblyInput
- Input data for Dirichlet boundary condition assembly.
-
template<class FES>class LinearFormAssemblyInput
- Input data for linear form assembly.
-
template<class ... Ts>class LinearFormTupleAssemblyInput
- Input data for tuple of linear forms assembly.
-
template<class U, class V>class MPI<Rodin::PETSc::Math::LinearSystem, Rodin::Variational::Problem<Rodin::PETSc::Math::LinearSystem, U, V>>
- MPI-parallel assembly of a single-variable PETSc problem.
-
template<class U1, class U2, class U3, class ... Us>class MPI<Rodin::PETSc::Math::LinearSystem, Rodin::Variational::Problem<Rodin::PETSc::Math::LinearSystem, U1, U2, U3, Us...>>
- MPI-parallel assembly of a multi-variable PETSc problem.
-
template<class Solution, class TrialFES, class TestFES>class MPI<::Mat, Variational::BilinearForm<Solution, TrialFES, TestFES, ::Mat>>
- MPI-parallel assembly of a PETSc matrix from a bilinear form.
-
template<class FES>class MPI<::Vec, Variational::LinearForm<FES, ::Vec>>
- MPI-parallel assembly of a PETSc vector from a linear form.
- class MPIIteration
- Iteration helper over a region of an MPI mesh shard.
-
template<class Solution, class TrialFES, class TestFES>class OpenMP<::Mat, Variational::BilinearForm<Solution, TrialFES, TestFES, ::Mat>>
- OpenMP-parallel assembly of a PETSc matrix from a bilinear form.
-
template<class FES>class OpenMP<::Vec, Variational::LinearForm<FES, ::Vec>>
- OpenMP-parallel assembly of a PETSc vector from a linear form.
-
template<class Solution, class TrialFES, class TestFES>class OpenMP<Math::SparseMatrix<typename FormLanguage::Dot<typename FormLanguage::Traits<TrialFES>::ScalarType, typename FormLanguage::Traits<TestFES>::ScalarType>::Type>, Variational::BilinearForm<Solution, TrialFES, TestFES, Math::SparseMatrix<typename FormLanguage::Dot<typename FormLanguage::Traits<TrialFES>::ScalarType, typename FormLanguage::Traits<TestFES>::ScalarType>::Type>>>
- OpenMP assembly of the Math::
SparseMatrix associated to a BilinearFormBase object. -
template<class FES>class OpenMP<Math::Vector<typename FormLanguage::Traits<FES>::ScalarType>, Variational::LinearForm<FES, Math::Vector<typename FormLanguage::Traits<FES>::ScalarType>>>
- OpenMP assembly of the Math::
Vector associated to a LinearForm object. -
template<class U, class V>class OpenMP<Rodin::PETSc::Math::LinearSystem, Rodin::Variational::Problem<Rodin::PETSc::Math::LinearSystem, U, V>>
- OpenMP-parallel assembly of a single-variable PETSc problem.
-
template<class U1, class U2, class U3, class ... Us>class OpenMP<Rodin::PETSc::Math::LinearSystem, Rodin::Variational::Problem<Rodin::PETSc::Math::LinearSystem, U1, U2, U3, Us...>>
- OpenMP-parallel assembly of a multi-variable PETSc problem.
-
template<class Solution, class TrialFES, class TestFES>class OpenMP<std::vector<Eigen::Triplet<typename FormLanguage::Dot<typename FormLanguage::Traits<TrialFES>::ScalarType, typename FormLanguage::Traits<TestFES>::ScalarType>::Type>>, Variational::BilinearForm<Solution, TrialFES, TestFES, std::vector<Eigen::Triplet<typename FormLanguage::Dot<typename FormLanguage::Traits<TrialFES>::ScalarType, typename FormLanguage::Traits<TestFES>::ScalarType>::Type>>>>
- OpenMP-based parallel assembly for bilinear forms.
-
template<>class OpenMPIteration<Geometry::Mesh<Context::Local>>
- OpenMP-based parallel mesh iteration for multi-threaded assembly.
-
template<class ProblemBody, class TrialFunction, class TestFunction>class ProblemAssemblyInput<ProblemBody, TrialFunction, TestFunction>
- Input data for complete problem assembly.
-
template<class LinearAlgebraType, class Operand>class Sequential
- Sequential (single-threaded) assembly implementation.
-
template<class U, class V>class Sequential<Rodin::PETSc::Math::LinearSystem, Rodin::Variational::Problem<Rodin::PETSc::Math::LinearSystem, U, V>>
- Sequential assembly of a single-variable PETSc problem.
-
template<class U1, class U2, class U3, class ... Us>class Sequential<Rodin::PETSc::Math::LinearSystem, Rodin::Variational::Problem<Rodin::PETSc::Math::LinearSystem, U1, U2, U3, Us...>>
- Sequential assembly of a multi-variable PETSc problem.
-
template<class Solution, class TrialFES, class TestFES>class Sequential<::Mat, Variational::BilinearForm<Solution, TrialFES, TestFES, ::Mat>>
- Sequential assembly of a PETSc matrix from a bilinear form.
-
template<class FES>class Sequential<::Vec, Variational::LinearForm<FES, ::Vec>>
- Sequential assembly of a PETSc vector from a linear form.
-
template<class Mesh>class SequentialIteration
- Sequential mesh iteration strategy.
-
template<>class SequentialIteration<Geometry::Mesh<Context::Local>>
- Sequential mesh iteration for single-threaded assembly.
Functions
-
template<class TrialFES, class TestFES>BilinearFormAssemblyInput(const TrialFES&, const TestFES&, FormLanguage::
List<Variational:: LocalBilinearFormIntegratorBase<decltype(std::declval<typename FormLanguage:: Traits<TrialFES>::ScalarType>()*std::declval<typename FormLanguage:: Traits<TestFES>::ScalarType>())>, &, FormLanguage:: List<Variational:: GlobalBilinearFormIntegratorBase<decltype(std::declval<typename FormLanguage:: Traits<TrialFES>::ScalarType>()*std::declval<typename FormLanguage:: Traits<TestFES>::ScalarType>())>, &) -> BilinearFormAssemblyInput< TrialFES, TestFES > - Template argument deduction guide for BilinearFormAssemblyInput.
-
OpenMPIteration(const Geometry::
Mesh<Context:: Local>& mesh, const Geometry:: Region&) -> OpenMPIteration< Geometry::Mesh< Context::Local > > - Template argument deduction guide for OpenMPIteration.
-
SequentialIteration(const Geometry::
Mesh<Context:: Local>& mesh, const Geometry:: Region&) -> SequentialIteration< Geometry::Mesh< Context::Local > > - Template argument deduction guide for SequentialIteration.
Function documentation
#include <Rodin/Assembly/Input.h>
template<class TrialFES, class TestFES>
Rodin:: Assembly:: BilinearFormAssemblyInput(const TrialFES&,
const TestFES&,
FormLanguage:: List<Variational:: LocalBilinearFormIntegratorBase<decltype(std::declval<typename FormLanguage:: Traits<TrialFES>::ScalarType>()*std::declval<typename FormLanguage:: Traits<TestFES>::ScalarType>())>,
&,
FormLanguage:: List<Variational:: GlobalBilinearFormIntegratorBase<decltype(std::declval<typename FormLanguage:: Traits<TrialFES>::ScalarType>()*std::declval<typename FormLanguage:: Traits<TestFES>::ScalarType>())>,
&) -> BilinearFormAssemblyInput< TrialFES, TestFES >
Template argument deduction guide for BilinearFormAssemblyInput.