Skip to content

Prior Builder

builder

Builders for prior objects from lower-level components.

Classes:

Name Description
BilaplacianPriorSettings

Settings for the bilaplacian prior builder.

BilaplacianPriorBuilder

Builder for a Bilaplacian prior.

ls_prior.builder.BilaplacianPriorSettings dataclass

Settings for the bilaplacian prior builder.

This dataclass collects all configuration options required to set up a biplaplacian prior using the BilaplacianPriorBuilder class. The builder distributes these settings to the respective components that are assembled within the builder.

Attributes:

Name Type Description
mesh dlx.mesh.Mesh

Dolfinx mesh on which the prior is defined.

mean_vector np.ndarray[tuple[int], np.dtype[np.float64]]

Mean vector of the prior, vertex-based representation.

kappa Real

Parameter \(\kappa\) in the SPDE formulation of the prior

tau Real

Parameter \(\tau\) in the SPDE formulation of the prior

robin_const Real | None

Robin boundary condition constant. If None, homogeneous Neumann boundary conditions are applied. Defaults to None.

seed int

Random seed for the internal random number generator. Defaults to 0.

fe_data tuple[str, int]

Finite element type and degree used for the function space setup. Defaults to ("CG", 1).

cg_relative_tolerance Real | None

Relative tolerance for the CG solver used in the application of the precision operator. If None, the default PETSc tolerance is used. Defaults to None.

cg_absolute_tolerance Real | None

Absolute tolerance for the CG solver used in the application of the precision operator. If None, the default PETSc tolerance is used. Defaults to None.

cg_max_iterations int | None

Maximum number of iterations for the CG solver used in the application of the precision operator. If None, the default PETSc value is used Defaults to None.

amg_relative_tolerance Real | None

Relative tolerance for the AMG solver used in the application of the covariance operator and its factorization. If None, the default PETSc tolerance is used. Defaults to None.

amg_absolute_tolerance Real | None

Absolute tolerance for the AMG solver used in the application of the covariance operator and its factorization. If None, the default PETSc tolerance is used. Defaults to None.

amg_max_iterations int | None

Maximum number of iterations for the AMG solver used in the application of the covariance operator and its factorization. If None, the default PETSc value is used. Defaults to None.

ls_prior.builder.BilaplacianPriorBuilder

Builder for a Bilaplacian prior.

Specific builder class for a bilaplacian prior, i.e. the distribution that arises from solving the SPDE \(\tau(\kappa^2 - \Delta)^2 m = W\).

Methods:

Name Description
build

Build the Bilaplacian prior.

__init__

__init__(settings: BilaplacianPriorSettings) -> None

Initialize the builder with the given settings.

Parameters:

Name Type Description Default
settings BilaplacianPriorSettings

Settings for the Bilaplacian prior.

required

build

build() -> prior.Prior

Build the Bilaplacian prior.

Internally, this method assembles all dolfinx structures, composes a hierarchy of PETScComponents, wraps them in InterfaceComponents and hands them to the Prior class.

Returns:

Type Description
prior.Prior

prior.Prior: The constructed Bilaplacian prior.

_build_fem_structures

_build_fem_structures() -> tuple[PETSc.Mat, PETSc.Mat, PETSc.Mat, PETSc.Mat, fem.FEMConverter]

Assemble FEM data structures, i.e. FEM Matrices and FEMConverter object.

Returns:

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

tuple[PETSc.Mat, PETSc.Mat, PETSc.Mat, fem.FEMConverter]: mass matrix \(M\), SPDE matrix \(A\), block-diagonal matrix \(\widehat{M}_w\), DoF map matrix \(L\), FEMConverter object.

_build_components

_build_components(
    mass_matrix: PETSc.Mat,
    spde_matrix: PETSc.Mat,
    block_diagonal_matrix: PETSc.Mat,
    dof_map_matrix: PETSc.Mat,
) -> tuple[components.PETScComponent, components.PETScComponent, components.PETScComponent]

Build PETSc components for bilaplacian prior from FEM matrices.

The main components for the prior are:

  1. Precision operator: \(\mathcal{C}^{-1} = A M^{-1} A\)
  2. Covariance operator: \(\mathcal{C} = A^{-1} M A^{-1}\)
  3. Sampling factor: \(\widehat{\mathcal{C}} = A^{-1} L^T \widehat{M}_e\)

Parameters:

Name Type Description Default
mass_matrix PETSc.Mat

Mass matrix \(M\).

required
spde_matrix PETSc.Mat

SPDE matrix \(A\).

required
block_diagonal_matrix PETSc.Mat

Block-diagonal matrix \(\widehat{M}_e\).

required
dof_map_matrix PETSc.Mat

DoF map matrix \(L\).

required

Returns:

Type Description
tuple[components.PETScComponent, components.PETScComponent, components.PETScComponent]

tuple[components.PETScComponent, components.PETScComponent, components.PETScComponent]: Precision operator \(\mathcal{C}^{-1}\), covariance operator \(\mathcal{C}\), sampling factor \(\widehat{\mathcal{C}}\).

_build_interfaces

_build_interfaces(
    precision_operator: components.PETScComponent,
    covariance_operator: components.PETScComponent,
    sampling_factor: components.PETScComponent,
) -> tuple[
    components.InterfaceComponent, components.InterfaceComponent, components.InterfaceComponent
]

Wrap PETSc components in interface components.

Parameters:

Name Type Description Default
precision_operator components.PETScComponent

Precision operator \(\mathcal{C}^{-1}\).

required
covariance_operator components.PETScComponent

SPDE matrix \(A\).

required
sampling_factor components.PETScComponent

Sampling factor \(\widehat{\mathcal{C}}\).

required

Returns:

Type Description
tuple[components.InterfaceComponent, components.InterfaceComponent, components.InterfaceComponent]

tuple[components.InterfaceComponent, components.InterfaceComponent, components.InterfaceComponent]: Wrapped interface components.