irksome.backends package¶
Submodules¶
irksome.backends.dolfinx module¶
DOLFINx backend for Irksome
- class irksome.backends.dolfinx.Constant(value)[source]¶
Bases:
ScalarValueInitialise.
- class irksome.backends.dolfinx.DirichletBC(*args: Any, **kwargs: Any)[source]¶
Bases:
DirichletBCCreate an Irksome compatible DirichletBC from an existing DOLFINx bc.
Note
This is not a user-facing class, but rather an internal class used for BC reconstruction in time-stepping problems. Use:
irksome.backends.dolfinx.dirichletbc.- Parameters:
g – The boundary condition expression
dofs – An array of degree-of-freedom indices in V
V – The space to construct the BC on.
- class irksome.backends.dolfinx.FloatConstantFunction(*args: Any, **kwargs: Any)[source]¶
Bases:
Function
- irksome.backends.dolfinx.Function(V: FunctionSpace | MixedFunctionSpace, name=None)[source]¶
Create a function in the backend language.
- class irksome.backends.dolfinx.LinearProblem(*args: Any, **kwargs: Any)[source]¶
Bases:
LinearProblemOverloaded linear problem that pack BCs before solving
- class irksome.backends.dolfinx.ListTensor(*expressions)[source]¶
Bases:
ListTensorA list tensor that exposes subfunctions for DOLFINx functions
Initialise.
- property subfunctions¶
- class irksome.backends.dolfinx.MeshConstant(msh)[source]¶
Bases:
object- Constant(val=0.0) Coefficient[source]¶
- class irksome.backends.dolfinx.NonlinearProblem(*args: Any, **kwargs: Any)[source]¶
Bases:
NonlinearProblemOverloaded nonlinear problem that pack BCs before solving.
Does each Newton iteration or line search step by overriding the SNES function and Jacobian assembly routines to pack BCs before assembly.
- property snes: SNES¶
- irksome.backends.dolfinx.TestFunction(function_space)[source]¶
Create a test function in the backend language.
This is required as
ufl.TestFunctionsreturns a tuple of test functions, which we need to convert to a tensor for consistency with the Firedrake backend.
- irksome.backends.dolfinx.TrialFunction(function_space)[source]¶
Create a trial function in the backend language.
This is required as
ufl.TrialFunctionsreturns a tuple of trial functions, which we need to convert to a tensor for consistency with the Firedrake backend.
- irksome.backends.dolfinx.assemble(expr: Expr | float)[source]¶
Assemble a UFL expression in the backend language.
- irksome.backends.dolfinx.create_bounds_constrained_bc(V, g, sub_domain, bounds, solver_parameters=None)[source]¶
- irksome.backends.dolfinx.create_variational_problem(F, u, bcs=None, aP=None, **kwargs)[source]¶
Create a variational problem.
- irksome.backends.dolfinx.create_variational_solver(problem: dolfinx.fem.petsc.LinearProblem | dolfinx.fem.petsc.NonlinearProblem, **kwargs)[source]¶
Create a variational solver that uses PETSc SNES or KSP.
- irksome.backends.dolfinx.dirichletbc(value: ufl.core.expr.Expr, dofs: npt.NDArray[np.int32], V: dolfinx.fem.FunctionSpace | None = None) DirichletBC[source]¶
Overloaded DirichletBC so that we can reconstruct BCs with UFL expressions.
Note
This class is user-facing.
- Parameters:
value – A UFL expression representing the boundary condition.
dofs – An array of degree-of-freedom indices in V where the BC should be applied.
V – The function space on which the BC applies. It can be a subspace of a mixed/blocked space.
- irksome.backends.dolfinx.extract_dtype(expr: Expr) DTypeLike[source]¶
Extract the dtype from an expression.
Looks for any constants or coefficients and returning their dtype. This is necessary for determining which DOLFINx DirichletBC constructor to use when packing UFL expressions into DOLFINx Expressions for use in BC reconstruction.
- irksome.backends.dolfinx.extract_function(expr) tuple[bool, dolfinx.fem.Function | None][source]¶
Recursively extract a Function from nested UFL expressions.
- Returns:
A tuple (is_real, func) where is_real=True means func is on RealElement (a constant)
Note
Returns (False, None) for RealElement functions so they get handled as scalars by extract_scalar_value, preserving any multipliers.
- irksome.backends.dolfinx.extract_linear_combination(expr: Expr) list[tuple[float, dolfinx.fem.Function]][source]¶
Extract (weight, function) pairs from a UFL linear combination.
Analyzes expressions like: 0.5*u + 0.3*v + 0.2*w Returns: [(0.5, u), (0.3, v), (0.2, w)]
- Parameters:
expr – UFL expression (Sum, Product, or single Function)
- Returns:
List of (weight, function) tuples
- irksome.backends.dolfinx.extract_scalar_value(scalar_expr)[source]¶
Extract float from a scalar UFL expression.
- irksome.backends.dolfinx.function_iadd(func: dolfinx.fem.Function, expr: Expr) None[source]¶
Add a UFL expression to a DOLFINx Function in-place (func += expr).
Extracts the linear combination structure and adds terms directly using PETSc vector operations. Works with mixed spaces and subfunction views.
- Parameters:
func – DOLFINx Function or subfunction view to update in-place
expr – UFL expression (linear combination of Functions)
function_iadd(u, 0.5 * v + 0.3 * w) # equivalent to u += 0.5*v + 0.3*w
- irksome.backends.dolfinx.function_iadd_method(self, other)[source]¶
In-place addition for DOLFINx functions (self += other).
This method is monkey-patched onto
dolfinx.fem.Functionto enable Firedrake-style arithmetic operations.- Parameters:
self – The function to be modified
other – The expression to add. Can be a UFL expression, Function, Constant, or scalar.
- Returns:
self (for chaining)
- irksome.backends.dolfinx.function_subfunctions(self)[source]¶
Get subfunctions for a DOLFINx function, which may be in a mixed space.
This function is called when updating u0 in time-stepping problems.
- irksome.backends.dolfinx.getNullspace(V, Vbig, num_stages, nullspace)[source]¶
Computes the nullspace for a multi-stage method.
- Parameters:
V – The
ufl.FunctionSpaceon which the original time-dependent PDE is posed.Vbig – The multi-stage
ufl.MixedFunctionSpacefor the stage problemnum_stages – The number of stages in the RK method
nullspace – The nullspace for the original problem.
On output, we produce a PETSc nullspace defining the nullspace for the multistage problem.
- irksome.backends.dolfinx.get_function_space(u: Coefficient | Argument) FunctionSpace[source]¶
- irksome.backends.dolfinx.get_mesh_constant(MC: MeshConstant | None) Expr[source]¶
- irksome.backends.dolfinx.get_stage_space(V: FunctionSpace, num_stages: int) FunctionSpace[source]¶
- irksome.backends.dolfinx.get_stages(V: dolfinx.fem.FunctionSpace, num_stages: int) ListTensor[source]¶
Given a function space for a single time-step, get a duplicate of this space, repeated num_stages times.
- Parameters:
V – Space for single step
num_stages – Number of stages
- Returns:
A coefficient in the new function space
- irksome.backends.dolfinx.invalidate_jacobian(solver: dolfinx.fem.petsc.LinearProblem)[source]¶
Invalidate the Jacobian matrix in the backend language.
irksome.backends.firedrake module¶
Firedrake backend for Irksome
- class irksome.backends.firedrake.BoundsConstrainedDirichletBC(V, g, sub_domain, bounds, solver_parameters=None)[source]¶
Bases:
DirichletBC- property function_arg¶
The value of this boundary condition.
- class irksome.backends.firedrake.MeshConstant(msh: Mesh)[source]¶
Bases:
object- Constant(val=0.0) Coefficient[source]¶
- irksome.backends.firedrake.create_bounds_constrained_bc(V, g, sub_domain, bounds, solver_parameters=None)[source]¶
- irksome.backends.firedrake.create_variational_problem(F, u, bcs=None, J=None, Jp=None, **kwargs)[source]¶
- irksome.backends.firedrake.extract_bcs(bcs: Any) tuple[Any][source]¶
Return an iterable of boundary conditions on the residual form
- irksome.backends.firedrake.getNullspace(V, Vbig, num_stages, nullspace)[source]¶
Computes the nullspace for a multi-stage method.
- Parameters:
V – The
firedrake.FunctionSpaceon which the original time-dependent PDE is posed.Vbig – The multi-stage
firedrake.FunctionSpacefor the stage problemnum_stages – The number of stages in the RK method
nullspace – The nullspace for the original problem.
On output, we produce a
firedrake.MixedVectorSpaceBasisdefining the nullspace for the multistage problem.
- irksome.backends.firedrake.get_function_space(u: Coefficient) FunctionSpace[source]¶
- irksome.backends.firedrake.get_mesh_constant(MC: MeshConstant | None)[source]¶
- irksome.backends.firedrake.get_stage_space(V: FunctionSpace, num_stages: int) FunctionSpace[source]¶
- irksome.backends.firedrake.get_stages(V: FunctionSpace, num_stages: int) Function[source]¶
Given a function space for a single time-step, get a duplicate of this space, repeated num_stages times.
- Parameters:
V – Space for single step
num_stages – Number of stages
- Returns:
A coefficient in the new function space