Skip to content

FEM-specific Functionalilty

fem

This module contains all FEM specific functionality in this package, based on dolfinx.

Classes:

Name Description
FEMConverter

Converter between vertex based data and DoF representation on a dolfinx function space.

FEMMatrixFactorizationAssembler

Assembler for the rectangular factorization of an FEM matrix.

Functions:

Name Description
generate_forms

Generate variational forms for the mass matrix and SPDE system matrix.

ls_prior.fem.ffi module-attribute

ffi = cffi.FFI()

ls_prior.fem.FEMConverter

Converter between vertex based data and DoF representation on a dolfinx function space.

This class connects the representation of arrays in dolfinx on a specified function space with a vertex-based viewpoint. Think of it as the adapter required for the prior to communicate with external components. The underlying idea is that such outside components only see the computational mesh, and define discrete data over the vertices of that mesh. An FEMConverter object takes such data structures and inerpolates them to the provided function space. On the other hand, it can interpolate any data defined on the DoFs of the underlying function space onto the vertices of the mesh. Internally, the FEMConverter assigns vertex based input data to the DoFs of a P1 function space (whose degrees of freedom are exactly the vertices). It subsequently utilizes dolfinx's efficient interpolation between function spaces. On the other hand, data from some function space is interpolated to a P1 space, and subsequently extracted to vertex values.

Methods:

Name Description
convert_vertex_values_to_dofs

Convert vertex based data to DoF representation

convert_dofs_to_vertex_values

Convert DoF based data to vertex representation

__init__

__init__(function_space: dlx.fem.FunctionSpace) -> None

Initialize the converter for a given dolfinx function space.

The function space implicitly carries the mesh, and thus the vertices we want to convert from/to.

Parameters:

Name Type Description Default
function_space dlx.fem.FunctionSpace

Function space in which degrees of freedom lie.

required

convert_vertex_values_to_dofs

convert_vertex_values_to_dofs(
    vertex_values: np.ndarray[tuple[int], np.dtype[np.float64]],
) -> np.ndarray[tuple[int], np.dtype[np.float64]]

Convert vertex based data to DoF representation.

Parameters:

Name Type Description Default
vertex_values np.ndarray[tuple[int], np.dtype[np.float64]]

Array of data defined on the vertices of the underlying computational mesh.

required

Raises:

Type Description
ValueError

Checks that vertex_values has the same dimension as a P1 function space on the underlying mesh.

Returns:

Type Description
np.ndarray[tuple[int], np.dtype[np.float64]]

np.ndarray[tuple[int], np.dtype[np.float64]]: Input data interpolated to the required function space DoFs, copied.

convert_dofs_to_vertex_values

convert_dofs_to_vertex_values(
    dof_values: np.ndarray[tuple[int], np.dtype[np.float64]],
) -> np.ndarray[tuple[int], np.dtype[np.float64]]

Convert DoF based data to vertex representation.

Parameters:

Name Type Description Default
dof_values np.ndarray[tuple[int], np.dtype[np.float64]]

Array of data defined on the DoFs of the underlying function space.

required

Raises:

Type Description
ValueError

Checks that the dimension of the input vector matches the number of DoFs on the underlying function space, copied.

Returns:

Type Description
np.ndarray[tuple[int], np.dtype[np.float64]]

np.ndarray[tuple[int], np.dtype[np.float64]]: Array of data interpolated to the vertices of the underlying mesh.

ls_prior.fem.FEMMatrixFactorizationAssembler

Assembler for the rectangular factorization of an FEM matrix.

This class provides the functionality for efficient factorization of a finite element matrix \(M\), i.e. it implements the (sparse representation of an) assembly of a rectangular matrix \(\widehat{M}\) s.th. \(M = \widehat{M}\widehat{M}^T\). The factorization exploits the characteristics of standard finite element assembly procedures. This assembly is typically done locally for the contribution of each mesh cell. The contributions are then gathered into a global matrix, with overlap at indices of vertices that are shared between cells. Now let \(N\) be the number of degrees of freedom of the finite element space, \(M\) the number of cells in the underlying mesh, and \(N_e\) the number of degrees of freedom per cell. The matrix \(M\) clearly has size \((N, N)\). Importantly, it can be decomposed as \(M=L^T M_e L\). Here, the Matrix \(M_e\in \mathbb{R}^{MN_e}\) is block-diagonal, containing the local FEM matrix contributions over all cells in each block. \(L\in \mathbb{R}^{MN_e \times N}\) is a sparse matrix that maps the local cell DoFs in each block \(M_e\) to their respective global DoFs. We can efficiently compute the Cholesky factor \(\widehat{M}_e\) of \(M_e\) by computing the cholesky factorization of each individual block. In particular, \(\widehat{M}_e\) is still sparse, as opposed to standard Cholesky factors of sparse matrices. \(\widehat{M}_e\) and \(L\) now represent a sparse factorization \(\widehat{M} = L^T \widehat{M}_e\) as above. The two matrices shouldn't be multiplied directly to avoid fill-in effects. However, matrix-vector products can be computed efficiently, by first multiplying with \(\widehat{M}_e\), and the result with \(L\). For further information on the factorization procedure, we refer to this publication.

Warning

The assembly procedure uses internals of dolfinx, which might be subject to change in future versions.

Methods:

Name Description
assemble

Assemble the block factorization.

__init__

__init__(mesh: dlx.mesh.Mesh, function_space: dlx.fem.FunctionSpace, form: ufl.Form) -> None

Initialize the block factorization assembler.

Parameters:

Name Type Description Default
mesh dlx.mesh.Mesh

Dolfinx Mesh object representing the computational domain.

required
function_space dlx.fem.FunctionSpace

Function space of the FEM problem.

required
form ufl.Form

Weak form of the FEM matrix to factorize.

required

assemble

assemble() -> tuple[PETSc.Mat, PETSc.Mat]

Assemble the rectangular matrix factorization.

Returns:

Type Description
tuple[PETSc.Mat, PETSc.Mat]

tuple[PETSc.Mat, PETSc.Mat]: PETSc Matrices \(\widehat{M}_e\) and \(L\) representing the rectangular factorization of an FEM matrix.

_init_assembly_kernel

_init_assembly_kernel(mpi_communicator: MPI.Comm, form: ufl.Form) -> cffi.FFI.CData

Initialize the dolfinx assembly kernel.

Kernel objects in dolfinx generate local matrix contributions over a single mesh cell, corresponding to a provided weak form. They are typically called internally during the assembly process of FEM system matrices. Under the hood, this method uses the ffcx jit compiler, and returns the resulting kernel's tabulate_tensor method.

Warning

This method uses internals of dolfinx, which might be subject to change in future versions.

Parameters:

Name Type Description Default
mpi_communicator MPI.Comm

MPI communicator for the kernel to use.

required
form ufl.Form

Weak form to compile per cell.

required

Returns:

Type Description
cffi.FFI.CData

cffi.FFI.CData: Compiled callable, returning the local FEM matrix contribution for a a single mesh cell.

_set_up_petsc_mats

_set_up_petsc_mats() -> tuple[PETSc.Mat, PETSc.Mat]

Initialize the PETSc matrices \(M_e\) and \(L\) for the assembly process.

Returns:

Type Description
tuple[PETSc.Mat, PETSc.Mat]

tuple[PETSc.Mat, PETSc.Mat]: Empty block diagonal matrix \(M_e\) and local-to-global DoF matrix \(L\).

_insert_in_block_diagonal_matrix

_insert_in_block_diagonal_matrix(
    global_ind: np.integer,
    cell_matrix: np.ndarray[tuple[int, int], np.dtype[np.float64]],
    block_diagonal_matrix: PETSc.Mat,
) -> None

Insert FEM matrix Cholesky factor of single mesh cell into block-diagonal matrix \(M_e\).

Insertion is done in-place.

Parameters:

Name Type Description Default
global_ind np.integer

Global index of the current mesh cell.

required
cell_matrix np.ndarray[tuple[int, int], np.dtype[np.float64]]

Local cell matrix contribution values.

required
block_diagonal_matrix PETSc.Mat

Global matrix \(M_e\) to insert into.

required

_insert_in_dof_map_matrix

_insert_in_dof_map_matrix(
    global_ind: np.integer,
    global_cell_dofs: np.ndarray[tuple[int], np.dtype[np.integer]],
    dof_map_matrix: PETSc.Mat,
) -> None

Insert global DoFs of a single mesh cell into local-to-global DoF matrix \(L\).

Insertion is done in-place.

Parameters:

Name Type Description Default
global_ind np.integer

Global index of the current mesh cell.

required
global_cell_dofs np.ndarray[tuple[int], np.dtype[np.integer]]

Global indices of the DoFs in the current cell.

required
dof_map_matrix PETSc.Mat

Local-to-global DoF matrix \(L\).

required

_assemble_matrices_over_cells

_assemble_matrices_over_cells(block_diagonal_matrix: PETSc.Mat, dof_map_matrix: PETSc.Mat) -> None

Assemble the matrices \(M_e\) and \(L\).

The method loops over all cells and invokes the kernel for the local FEM matrix contributions. It further computes the Cholesky factor of each local matrix, and inserts the result and global index mapping into the matrices \(M_e\) and \(L\), respectively.

Parameters:

Name Type Description Default
block_diagonal_matrix PETSc.Mat

Block diagonal matrix \(M_e\) for cell Cholesky factors.

required
dof_map_matrix PETSc.Mat

Local-to-global DoF matrix \(L\).

required

ls_prior.fem.generate_forms

generate_forms(
    function_space: dlx.fem.FunctionSpace,
    kappa: Annotated[Real, Is[lambda x: x > 0]],
    tau: Annotated[Real, Is[lambda x: x > 0]],
    robin_const: Real | None = None,
) -> tuple[ufl.Form, ufl.Form]

Construct dolfinx forms for the mass matrix and SPDE system matrix.

This method constructs dolfinx variational forms that resemble the mass matrix contribution and SPDE system matrix, i.e. the left-hand-side of the SPDE generating the random field. More specifically, let \(\Omega\) be the domain of interest, and \(\phi\) both the trial and test function (defined over the same function space). The mass matrix contribution is then given as \((\phi, \phi)_{L^2(\Omega)}\), and the total SPDE matrix contribution is

\[ \begin{equation*} \kappa^2 \tau (\phi, \phi)_{L^2(\Omega)} + \tau (\nabla \phi, \nabla \phi)_{L^2(\Omega)} + \beta (\phi, \phi)_{L^2(\partial\Omega)} \end{equation*} \]

\(\kappa\) and \(\tau\) correspond to the parameters kappa and tau for the parameterization of the field. \(\beta\) is the optional robin_const parameter that enforces Robin boundary conditions of the form \(\nabla \phi \cdot n + \beta \phi = 0\), instead of homogeneous Neumann boundary conditions.

Parameters:

Name Type Description Default
function_space dlx.fem.FunctionSpace

dolfinx function space to construct forms over.

required
kappa Real

Parameter \(\kappa\) of the prior field.

required
tau Real

Parameter \(\tau\) of the prior field.

required
robin_const Real | None

Parameter \(\beta\) of the prior field enforcing Robin boundary conditions. Defaults to None.

None

Returns:

Type Description
tuple[ufl.Form, ufl.Form]

tuple[ufl.Form, ufl.Form]: dolfinx forms for the mass matrix and SPDE system matrix contributions.