API reference#

class fenicsx_ii.Average(u: Argument, restriction_operator, restriction_space)[source]#

Averaging operator.

Parameters:
  • u – Argument to be averaged.

  • restriction_operator – Operator that restricts from the parent space to the averaged space.

  • restriction_space – Function space on the new mesh.

class fenicsx_ii.Circle(mesh: Mesh, radius: float | Callable[[ndarray[tuple[int, ...], dtype[floating]]], ndarray[tuple[int, ...], dtype[floating]]], degree: int = 15)[source]#

The average on the perimiter of a circle with a given radius.

Parameters:
  • mesh – The 1D mesh on which the circles are centered.

  • radius – The radius of the circles. Either a float or a function that takes in a list of points (shape (3, num_points)) and returns a list of radii (shape (num_points,)).

  • degree – The degree of the quadrature rule on the circle. Must be at least 3 for a sensible approximation.

quadrature(x0: ndarray[tuple[int, ...], dtype[floating]], normals: ndarray[tuple[int, ...], dtype[floating]]) tuple[ndarray[tuple[int, ...], dtype[floating]], ndarray[tuple[int, ...], dtype[floating]], ndarray[tuple[int, ...], dtype[floating]]][source]#

Compute quadrature points in physical space for circles with centers x0 and normals.

Parameters:
  • x0 – List of center points of the circles, shape (num_points, 3)

  • normals – List of normal vectors to the circles, shape (num_points, 3)

Returns:

Quadrature points and weights and scales for each weight in 3D.

class fenicsx_ii.Disk(mesh: Mesh, radius: float | Callable[[ndarray[tuple[int, ...], dtype[floating]]], ndarray[tuple[int, ...], dtype[floating]]], degree: int = 7)[source]#

The average on a disk with a given radius. The quadrature rule is based on Bojanov and Petrova (1998), https://doi.org/10.1007/s002110050358.

Parameters:
  • mesh – The 1D mesh on which the disks are centered.

  • radius – The radius of the disks. Either a float or a function that takes in a list of points (shape (3, num_points)) and returns a list of radii (shape (num_points,)).

  • degree – The degree of the quadrature rule on the disk.

quadrature(x0: ndarray[tuple[int, ...], dtype[floating]], normals: ndarray[tuple[int, ...], dtype[floating]]) tuple[ndarray[tuple[int, ...], dtype[floating]], ndarray[tuple[int, ...], dtype[floating]], ndarray[tuple[int, ...], dtype[floating]]][source]#

Compute quadrature points in physical space for circles with centers x0 and normals.

Parameters:
  • x0 – List of center points of the circles, shape (num_points, 3)

  • normals – List of normal vectors to the circles, shape (num_points, 3)

Returns:

Quadrature points and weights and scales for each weight in 3D.

class fenicsx_ii.LinearProblem(a: Form, L: Form, *, petsc_options_prefix: str, bcs: Sequence[DirichletBC] | None = None, u: Sequence[Function] | None = None, P: Form | None = None, petsc_options: dict | None = None, form_compiler_options: dict | None = None, jit_options: dict | None = None, entity_maps: Sequence[EntityMap] | None = None)[source]#

High-level class for solving a linear variational problem using a PETSc KSP.

Solves problems of the form \(a_{ij}(u, v) = f_i(v), i,j=0,\ldots,N\ \forall v \in V\) where \(u=(u_0,\ldots,u_N), v=(v_0,\ldots,v_N)\) using PETSc KSP as the linear solver.

Note

This high-level class automatically handles PETSc memory management. The user does not need to manually call .destroy() on returned PETSc objects.

property A: Mat#

Left-hand side matrix.

property P_mat: Mat#

Preconditioner matrix.

property b: Vec#

Right-hand side vector.

solve() Function | Sequence[Function][source]#

Solve the problem.

This method updates the solution u function(s) stored in the problem instance.

Note

The user is responsible for asserting convergence of the KSP solver e.g. problem.solver.getConvergedReason() > 0. Alternatively, pass "ksp_error_if_not_converged" : True in petsc_options to raise a PETScError on failure.

Returns:

The solution function(s).

property solver: KSP#

The PETSc KSP solver.

property u: Function | Sequence[Function]#

Solution function(s).

Note

The function(s) do not share memory with the solution vector x.

property x: Vec#

Solution vector.

Note

The vector does not share memory with the solution function(s) u.

class fenicsx_ii.PointwiseTrace(mesh: Mesh)[source]#

Trace of the 3D function on the line segment.

Parameters:

mesh – The 1D mesh on which the trace is taken.

class fenicsx_ii.Quadrature(name: str, points: numpy.ndarray[tuple[int, ...], numpy.dtype[numpy.floating]], weights: numpy.ndarray[tuple[int, ...], numpy.dtype[numpy.floating]], scales: numpy.ndarray[tuple[int, ...], numpy.dtype[numpy.floating]])[source]#
name: str#

Name of the quadrature rule

points: ndarray[tuple[int, ...], dtype[floating]]#

Points in the rule (physical space)

scales: ndarray[tuple[int, ...], dtype[floating]]#

Scaling factors for the quadrature points

weights: ndarray[tuple[int, ...], dtype[floating]]#

Weights of the quadrature points

class fenicsx_ii.ReductionOperator(*args, **kwargs)[source]#

Protocol for reduction operators from 3D to 1D.

compute_quadrature(cells: ndarray[tuple[int, ...], dtype[int32]], reference_points: ndarray[tuple[int, ...], dtype[floating]]) Quadrature[source]#

Function that should compute the Quadrature on the given cells and reference points.

Parameters:
  • cells – List of cell indices on which to compute the quadrature (local to process).

  • reference_points – List of reference points in the cells, shape (num_points, tdim).

Returns:

The quadrature rule on the given cells and reference points.

property num_points: int#

Number of quadrature points per cell.

fenicsx_ii.assemble_matrix(a: Form, bcs: Sequence[DirichletBC] | None = None, form_compiler_options: dict | None = None, jit_options: dict | None = None, entity_maps: Sequence[EntityMap] | None = None, A: Mat | None = None) Mat[source]#

Assemble a bi-linear ufl.Form into a petsc4py.PETSc.Mat.

The form might be a block form, in which case a nested matrix is returned.

Parameters:
  • a – Bi-linear UFL form to assemble.

  • bcs – List of Dirichlet boundary conditions to apply.

  • form_compiler_options – Options to pass to the form compiler.

  • jit_options – Options to pass to the JIT compiler.

  • A – An optional PETSc matrix to which the resulting matrix is added. If None, a new matrix is created.

fenicsx_ii.assemble_scalar(form: ~ufl.form.Form, form_compiler_options: dict | None = None, jit_options: dict | None = None, entity_maps: ~typing.Sequence[~dolfinx.mesh.EntityMap] | None = None, op: ~mpi4py.MPI.Op = <mpi4py.MPI.Op object>) inexact | float | complex[source]#

Assemble a rank-0 form into a scalar value. Perform necessary MPI communication.

Parameters:
  • form – The rank-0 form to assemble.

  • form_compiler_options – Options for the form compiler.

  • jit_options – Options for just-in-time compilation.

  • entity_maps – Entity maps for the mesh.

  • op – The MPI operation to use for communication (default is mpi4py.MPI.SUM).

Returns

The assembled scalar value.

fenicsx_ii.assemble_vector(L: Form, bcs: Sequence[DirichletBC] | None = None, form_compiler_options: dict | None = None, jit_options: dict | None = None, entity_maps: Sequence[EntityMap] | None = None, b: Vec | None = None) Vec[source]#

Assemble a ufl.Form into a petsc4py.PETSc.Vec.

The form might contain fenicsx_ii.Average operators, in which case the appropriate restriction operators are applied to the test and trial functions before creating the vector. The form might also be a block form, in which case a nested vector is created.

Parameters:
  • L – Linear UFL form to assemble.

  • bcs – List of Dirichlet boundary conditions to apply.

  • form_compiler_options – Options to pass to the form compiler.

  • jit_options – Options to pass to the JIT compiler.

  • entity_maps – Entity maps to use for assembly.

  • b – An optional PETSc vector to which the assembled vector is added. If None, a new vector is created.

fenicsx_ii.average_coefficients(form: Form)[source]#

Perform the average operation on all coefficients of a given form.

Note

This updates the underlying vector of each coefficient.

Parameters:

form – The form whose coefficients are averaged.

fenicsx_ii.create_interpolation_matrix(V: dolfinx.fem.FunctionSpace, K: dolfinx.fem.FunctionSpace, red_op: ReductionOperator, tol: float = 1e-08, use_petsc: bool = False) tuple['PETSc.Mat' | dolfinx.la.MatrixCSR, _im, _im][source]#

Create an interpolation matrix from V to K with a specific reduction operator applied to the interpolation points of K.

Parameters:
  • V – The function space to interpolate from.

  • K – The function space to interpolate to.

  • red_op – Reduction operator for each interpolation point on K.

  • tol – Tolerance for determining point ownership across processes.

  • use_petsc – Flag to indicate whether to use PETSc for the matrix.

Returns:

The interpolation matrix, as a DOLFINx built in matrix or a PETSc matrix. The accumulation of contributions from all processes has been accounted before returning the matrix.

fenicsx_ii.create_matrix(a: Form, form_compiler_options: dict | None = None, jit_options: dict | None = None, entity_maps: Sequence[EntityMap] | None = None) Mat[source]#

Create a petsc4py.PETSc.Mat from a ufl.Form.

The form might contain fenicsx_ii.Average operators, in which case the appropriate restriction operators are applied to the test and trial functions before creating the matrix. The form might also be a block form, in which case a nested matrix is created.

Parameters:
  • a – Bi-linear UFL form to assemble.

  • form_compiler_options – Options to pass to the form compiler.

  • jit_options – Options to pass to the JIT compiler.

fenicsx_ii.create_vector(L: Form, form_compiler_options: dict | None = None, jit_options: dict | None = None, entity_maps: Sequence[EntityMap] | None = None) Vec[source]#

Create a petsc4py.PETSc.Vec from a ufl.Form.

The form might contain fenicsx_ii.Average operators, in which case the appropriate restriction operators are applied to the test functions before creating the vector. The form might also be a block form, in which case a nested vector is created.

Parameters:
  • b – Linear UFL form to assemble.

  • form_compiler_options – Options to pass to the form compiler.

  • jit_options – Options to pass to the JIT compiler.