Laplacian Error Models#
- class DolfinLaplacianErrorModelConfigs(*, debug=False, verbose=False, output_dir='./_PYOED_RESULTS_', design=True, model=None, target='parameter', gamma=1.0, delta=0.5, mean=0.0, random_seed=None)[source]#
Bases:
ErrorModelConfigsConfigurations for
DolfinLaplacianErrorModel- Parameters:
verbose (bool) – a boolean flag to control verbosity of the object.
debug (bool) – a boolean flag that enables adding extra functionality in a debug mode
output_dir (str | Path) – the base directory where the output files will be saved.
design (None | bool | Sequence[bool] | ndarray) –
an experimental design to define active/inactive entries of the random variable (mean, variance/covariance matrix).
Note
If the design is None, it is set to all ones; that is everything is observed (default)
If the design is a binary vector ( or int dtype attributes with 0/1 entries) the mean, the covariance, and all random vectors are projected onto the space identified by the 1/True entries.
model (SimulationModel | None) – the simulation model from which the source grid is extracted
target (str) – either “parameter” or “state”. This defines whether to use model.parameter_dof or model.state_dof to create the finite element space to which the distribution is associated
gamma (float) – positive number (float) that controls the variance/uncertainty level
delta (float) – positive number (float) such that gamma / delta controls the correlation length (lengthscale),
mean (float | ndarray) – float or iterable of floats (numpy array) of the error model;
random_seed (int | None) – random seed used when the Gaussian model is initiated. This is useful for reproductivity. If None, random seed follows numpy.random.seed rules
- model: SimulationModel | None#
- __init__(*, debug=False, verbose=False, output_dir='./_PYOED_RESULTS_', design=True, model=None, target='parameter', gamma=1.0, delta=0.5, mean=0.0, random_seed=None)#
- class DolfinLaplacianErrorModel(configs=None)[source]#
Bases:
ErrorModel,RandomNumberGenerationMixinA class implementing the simplest Error model with a Laplacian-based covariance model. The covariance matrix takes the form \(C = P^{-1} M P^{-1}\), where \(P = \delta M + \gamma K\) is the precision (stiffness + mass), \(M\) is the FEM mass matrix, \(\gamma / \delta\) controls the correlation length (lengthscale), \(\gamma\) controls the variance/uncertainty level, and \(K\) is the stiffness matrix (discretised \(-\Delta\)).
Assembled using FEniCSx (dolfinx 0.10) / PETSc.
- Parameters:
configs (dict | DolfinLaplacianErrorModelConfigs | None) – settings of the error model.
- Raises:
ImportError: if dolfinx or ufl are not available/installed
TypeError: if mean type is not supported (scalars and iterables are supported)
ValueError: if no valid FEM space (DOF) Vh is found in the configurations dictionary
ValueError: if gamma or delta are non-positive
- validate_configurations(configs, raise_for_invalid=False)[source]#
Validate the passed configurations for the DolfinLaplacian error model
- Parameters:
configs (dict | DolfinLaplacianErrorModelConfigs) – configurations to validate.
raise_for_invalid (bool) – raise an exception if the configurations are invalid
- Returns:
True if the configurations are valid, otherwise False
- Raises:
PyOEDConfigsValidationError – if the configurations are invalid and raise_for_invalid is set to True.
- update_covariance_matrix(gamma, delta, validate=False)[source]#
Assemble the precision matrix \(P = \delta M + \gamma K\) and its KSP solver, then compute the dense Cholesky factor of M for sqrt operations.
- Parameters:
gamma – positive number controlling variance/uncertainty level
delta – positive number such that gamma / delta controls correlation length
validate (bool) – validate the covariance parameters before updating.
- update_mean(mean, validate=False)[source]#
Update the mean of the distribution.
- Parameters:
mean – mean of the distribution (scalar, numpy array, or PETSc.Vec)
validate (bool) – validate the mean before updating.
- create_noise_vector(return_np=True)[source]#
Create a vector compatible with the domain of the covariance sqrt operator. Since noise_size == state_size in the dolfinx implementation, this is the same as the state space.
- create_vector(return_np=True)[source]#
Create a vector compatible with the range of the underlying covariance operator.
- generate_white_noise(truncate=False, return_np=False)[source]#
Generate a standard normal random vector of size noise_size.
- add_noise(x, in_place=False)[source]#
Add random noise to the passed vector x.
- Parameters:
x – numpy array or PETSc.Vec of size equal to self.size
in_place (bool) – overwrite x if True.
- generate_noise(return_np=True)[source]#
Generate a zero-mean noise vector drawn from \(\mathcal{N}(0, C)\).
- covariance_sqrt_matvec(x, lower=True, return_np=True, in_place=False)[source]#
Apply the square root (or its transpose) of the covariance operator to x.
The covariance is \(C = P^{-1} M P^{-1}\), with square root \(C^{1/2} = P^{-1} L_M\) where \(M = L_M L_M^T\) (dense Cholesky).
- Parameters:
- Returns:
output vector of size state_size (lower=True) or noise_size (lower=False).
- covariance_matvec(x, return_np=True, in_place=False)[source]#
Apply the covariance \(C = P^{-1} M P^{-1}\) to x.
- Parameters:
x – input vector (numpy array or PETSc.Vec) of size self.size.
return_np (bool) – return numpy array (True) or PETSc.Vec (False).
- Returns:
output vector of size self.size.
- covariance_diagonal(method='exact', verbose=False)[source]#
Return a copy of the diagonal of the covariance matrix (the variances).
- Returns:
1D numpy array of size self.size.
- variance()[source]#
Same as
covariance_diagonal()
- covariance_inv_matvec(x, in_place=False, return_np=True)[source]#
Apply the precision (inverse covariance) to x via GMRES.
- Parameters:
x – input vector (numpy array or PETSc.Vec) of size self.size.
- Returns:
output vector of size self.size.
- covariance_trace(method='exact', sample_size=50, optimize_sample_size=False, min_sample_size=10, max_sample_size=100, distribution='Rademacher', accuracy=0.01)[source]#
Evaluate/Approximate the trace of the covariance matrix.
- covariance_logdet(method='exact', sample_size=50, optimize_sample_size=False, min_sample_size=10, max_sample_size=100, distribution='Rademacher', accuracy=0.01)[source]#
Evaluate/approximate the log-determinant of the covariance matrix.
- property model#
Return a reference to the underlying model
- property Vh#
Get the finite element function space
- property mean#
Get the distribution mean as a numpy array.
- property size#
Dimension of the underlying probability distribution
- property noise_size#
Dimension of the noise space (equals state_size in this implementation)
- property design#
Get a COPY of the design vector
- property active#
Flags (boolean 1d array-like) corresponding to active/inactive dimensions according to the experimental design. This vector is of exactly the same size as the observation vector (ignoring the experimental design).
Warning
Because the return type is bool, this vector can be used to slice the vector. To avoid any confustion, however, we decided to add a more clear property active_indexes which returns the indexes of the variable/vector corresponding to non-zero design values. Thus, one should always resort to using active_indexes in case extracting active entries in the variable (according to the design) is required.
- property active_indexes#
An array of the indexes of the observation vector that is active (associated with non-zero design) according to the experimental design.
- property active_size#
Dimension of the probability distribution projected onto the design space (only active/nonzero entries of the design)
- class DolfinBiLaplacianErrorModelConfigs(*, debug=False, verbose=False, output_dir='./_PYOED_RESULTS_', design=True, model=None, target='parameter', gamma=1.0, delta=0.5, mean=0.0, random_seed=None)[source]#
Bases:
DolfinLaplacianErrorModelConfigsConfigurations for
DolfinBiLaplacianErrorModel. Same asDolfinLaplacianErrorModelConfigs.- Parameters:
verbose (bool) – a boolean flag to control verbosity of the object.
debug (bool) – a boolean flag that enables adding extra functionality in a debug mode
output_dir (str | Path) – the base directory where the output files will be saved.
design (None | bool | Sequence[bool] | ndarray) –
an experimental design to define active/inactive entries of the random variable (mean, variance/covariance matrix).
Note
If the design is None, it is set to all ones; that is everything is observed (default)
If the design is a binary vector ( or int dtype attributes with 0/1 entries) the mean, the covariance, and all random vectors are projected onto the space identified by the 1/True entries.
model (SimulationModel | None) – the simulation model from which the source grid is extracted
target (str) – either “parameter” or “state”. This defines whether to use model.parameter_dof or model.state_dof to create the finite element space to which the distribution is associated
gamma (float) – positive number (float) that controls the variance/uncertainty level
delta (float) – positive number (float) such that gamma / delta controls the correlation length (lengthscale),
mean (float | ndarray) – float or iterable of floats (numpy array) of the error model;
random_seed (int | None) – random seed used when the Gaussian model is initiated. This is useful for reproductivity. If None, random seed follows numpy.random.seed rules
- __init__(*, debug=False, verbose=False, output_dir='./_PYOED_RESULTS_', design=True, model=None, target='parameter', gamma=1.0, delta=0.5, mean=0.0, random_seed=None)#
- class DolfinBiLaplacianErrorModel(configs=None)[source]#
Bases:
DolfinLaplacianErrorModelA class implementing a simple Error model with a Bi-Laplacian-based covariance model. The covariance matrix takes the form: \(\mathbf{C} = \mathbf{B}^{-1} \mathbf{M} \mathbf{B}^{-1}\), where \(\mathbf{B} = \delta \mathbf{M} + \gamma \mathbf{K}\) is the precision square root (FEM stiffness + mass), \(\gamma / \delta\) controls the correlation length (lengthscale), and \(\delta \gamma\) controls the variance/uncertainty level.
Assembled using FEniCSx (dolfinx 0.10) / PETSc.
- Parameters:
configs (DolfinBiLaplacianErrorModelConfigs | dict | None) – settings of the error model
- Raises:
ImportError: if dolfinx or ufl are not available/installed
TypeError: if mean or variance types are not supported
ValueError: if no valid FEM space (DOF) Vh is found in the configurations
ValueError: if gamma or delta are non-positive
- validate_configurations(configs, raise_for_invalid=False)[source]#
Validate the passed configurations.
- update_covariance_matrix(gamma, delta, validate=False)[source]#
Assemble the BiLaplacian precision square-root matrix \(B = \delta M + \gamma K\) and its KSP solver, then compute the dense Cholesky factor of M for sqrt operations.
- Parameters:
gamma – positive number (float) controlling variance/uncertainty level
delta – positive number (float) controlling correlation length
validate (bool) – validate the covariance parameters before updating.
- create_noise_vector(return_np=True)[source]#
Create a vector compatible with the noise domain (= state space).
- covariance_sqrt_matvec(x, lower=True, return_np=True, in_place=False)[source]#
Apply \(C^{1/2} = B^{-1} L_M\) (lower=True) or \((C^{1/2})^T = L_M^T B^{-1}\) (lower=False) to x, where \(M = L_M L_M^T\) (dense Cholesky).
- Parameters:
x – input vector (numpy array or PETSc.Vec).
lower (bool) – apply forward (True) or adjoint (False) sqrt.
- Returns:
output vector.
- create_Laplacian_error_model(model, target='parameter', gamma=1.0, delta=0.5, mean=0.0, random_seed=None)[source]#
A simple interface to instantiate an object of the
DolfinLaplacianErrorModelclass.- Parameters:
model – the simulation model from which the source grid is extracted
target – either “parameter” or “state”
gamma – positive number controlling variance/uncertainty level
delta – positive number controlling correlation length
mean – mean of the distribution (scalar or array)
random_seed – random seed for reproducibility
- Returns:
an instance of
DolfinLaplacianErrorModel
- create_BiLaplacian_error_model(model, target='parameter', gamma=1.0, delta=0.5, mean=0.0, random_seed=None)[source]#
A simple interface to instantiate an object of the
DolfinBiLaplacianErrorModelclass.- Parameters:
model – the simulation model from which the source grid is extracted
target – either “parameter” or “state”
gamma – positive number controlling variance/uncertainty level
delta – positive number controlling correlation length
mean – mean of the distribution (scalar or array)
random_seed – random seed for reproducibility
- Returns:
an instance of
DolfinBiLaplacianErrorModel