FEniCSPDE#
- class cuqipy_fenics.pde.FEniCSPDE(PDE_form, mesh, solution_function_space, parameter_function_space, dirichlet_bcs, adjoint_dirichlet_bcs=None, observation_operator=None, reuse_assembled=False, linalg_solve=None, linalg_solve_kwargs=None)#
Base class that represents PDEs problem defined in FEniCS. This class is not meant to be used directly, but rather it defines the API that should be implemented by subclasses for specific types of PDEs.
This class along with the cuqipy_fenics.geometry, interfaces FEniCS PDE models with CUQIpy library.
- Parameters:
PDE_form (callable or tuple of two callables) –
If passed as a callable: the callable returns the weak form of the PDE. The callable should take three arguments, the first argument is the parameter (input of the forward model), the second argument is the state variable (solution variable), and the third argument is the adjoint variable (the test variable in the weak formulation).
If passed as a tuple of two callables: the first callable returns the weak form of the PDE left hand side, and the second callable returns the weak form of the PDE right hand side. The left hand side callable takes the same three arguments as described above. The right hand side callable takes only the parameter and the adjoint variable as arguments. See the example below.
mesh (FEniCS mesh) – FEniCS mesh object that defines the discretization of the domain.
solution_function_space (FEniCS function space) – FEniCS function space object that defines the function space of the state variable (solution variable).
parameter_function_space (FEniCS function space) – FEniCS function space object that defines the function space of the Bayesian parameter (input of the forward model).
dirichlet_bcs (FEniCS Dirichlet boundary condition object or a list of them) – FEniCS Dirichlet boundary condition object(s) that define the Dirichlet boundary conditions of the PDE.
adjoint_dirichlet_bcs (FEniCS Dirichlet boundary condition object or a list of them, required if the gradient is to be computed) – FEniCS Dirichlet boundary condition object(s) that define the Dirichlet boundary conditions of the adjoint PDE.
observation_operator (python function handle, optional) –
Function handle of a python function that returns the observed quantity from the PDE solution. If not provided, the identity operator is assumed (i.e. the entire solution is observed). This python function takes as input the Bayesian parameter (input of the forward model) and the state variable (solution variable) as first and second inputs, respectively.
The returned observed quantity can be a ufl.algebra.Operator, FEniCS Function, np.ndarray, int, or float.
Examples of observation_operator are lambda m, u: m*u and lambda m, u: m*dl.sqrt(dl.inner(dl.grad(u),dl.grad(u)))
reuse_assembled (bool, optional) – Flag to indicate whether the assembled (and possibly factored) differential operator should be reused when the parameter is not changed. If True, the assembled matrices are reused. If False, the assembled matrices are not reused. Default is False.
linalg_solve (FEniCS linear solver object, optional)
linalg_solve_kwargs (dict, optional) – Dictionary of keyword arguments to be passed to the linear solver object.
Example
# Define mesh mesh = dl.UnitSquareMesh(20, 20) # Set up function spaces solution_function_space = dl.FunctionSpace(mesh, 'Lagrange', 2) parameter_function_space = dl.FunctionSpace(mesh, 'Lagrange', 1) # Set up Dirichlet boundaries def u_boundary(x, on_boundary): return on_boundary dirichlet_bc_expr = dl.Expression("0", degree=1) dirichlet_bcs = dl.DirichletBC(solution_function_space, dirichlet_bc_expr, u_boundary) # Set up PDE variational form def lhs_form(m,u,p): return ufl.exp(m)*ufl.inner(ufl.grad(u), ufl.grad(p))*ufl.dx def rhs_form(m,p): return - dl.Constant(1)*p*ufl.dx # Create the PDE object PDE = cuqipy_fenics.pde.SteadyStateLinearFEniCSPDE( (lhs_form, rhs_form), mesh, parameter_function_space=parameter_function_space, solution_function_space=solution_function_space, dirichlet_bcs=dirichlet_bcs)
- __init__(PDE_form, mesh, solution_function_space, parameter_function_space, dirichlet_bcs, adjoint_dirichlet_bcs=None, observation_operator=None, reuse_assembled=False, linalg_solve=None, linalg_solve_kwargs=None)#
Methods
__init__
(PDE_form, mesh, ...[, ...])assemble
(parameter)Assemble the PDE weak form
Compute gradient of the PDE weak form w.r.t.
observe
(PDE_solution)Apply observation operator on the PDE solution
solve
()Solve the PDE
Attributes
Get the PDE form
Get the forward solution of the PDE
Get the lhs form
Get the observation operator
Get the parameter of the PDE
Get the reuse_assembled flag
Get the rhs form