dolfinx.fem

Tools for assembling and manipulating finite element forms.

Functions

create_sparsity_pattern(a)

Create a sparsity pattern from a bilinear form.

class dolfinx.fem.Constant(domain, c: Union[ndarray, Sequence, float, complex])[source]

Bases: Constant

A constant with respect to a domain.

Parameters:
  • domain – DOLFINx or UFL mesh

  • c – Value of the constant.

property dtype: dtype
property value

The value of the constant

class dolfinx.fem.DirichletBC(bc)[source]

Bases: object

Representation of Dirichlet boundary condition which is imposed on a linear system.

Note

Dirichlet boundary conditions should normally be constructed using fem.dirichletbc() and not using this class initialiser. This class is combined with different base classes that depend on the scalar type of the boundary condition.

Parameters:
  • value – Lifted boundary values function.

  • dofs – Local indices of degrees of freedom in function space to which boundary condition applies. Expects array of size (number of dofs, 2) if function space of the problem, V is passed. Otherwise assumes function space of the problem is the same of function space of boundary values function.

  • V – Function space of a problem to which boundary conditions are applied.

property function_space: FunctionSpaceBase

The function space on which the boundary condition is defined

property g: Union[Function, Constant, np.ndarray]

The boundary condition value(s)

class dolfinx.fem.DofMap(dofmap: DofMap)[source]

Bases: object

Degree-of-freedom map

This class handles the mapping of degrees of freedom. It builds a dof map based on a ufcx_dofmap on a specific mesh.

property bs

Returns the block size of the dofmap

cell_dofs(cell_index: int)[source]

Cell local-global dof map

Parameters:

cell – The cell index

Returns:

Local-global dof map for the cell (using process-local indices)

property dof_layout

Layout of dofs on an element

property index_map

Index map that described the parallel distribution of the dofmap

property index_map_bs

Block size of the index map

property list

Adjacency list with dof indices for each cell

class dolfinx.fem.ElementMetaData(family: str, degree: int, shape: Optional[Tuple[int, ...]] = None, symmetry: Optional[bool] = None)[source]

Bases: NamedTuple

Data for representing a finite element

Parameters:
  • family – Element type.

  • degree – Polynomial degree of the element.

  • shape – Shape for vector/tensor valued elements that are constructed from blocked scalar elements (e.g., Lagrange).

  • symmetry – Symmetry option for blocked tensor elements.

Create new instance of ElementMetaData(family, degree, shape, symmetry)

degree: int

Alias for field number 1

family: str

Alias for field number 0

shape: Optional[Tuple[int, ...]]

Alias for field number 2

symmetry: Optional[bool]

Alias for field number 3

class dolfinx.fem.Expression(e: ufl.core.expr.Expr, X: np.ndarray, comm: Optional[_MPI.Comm] = None, form_compiler_options: Optional[dict] = None, jit_options: Optional[dict] = None, dtype: Optional[npt.DTypeLike] = None)[source]

Bases: object

Create DOLFINx Expression.

Represents a mathematical expression evaluated at a pre-defined set of points on the reference cell. This class closely follows the concept of a UFC Expression.

This functionality can be used to evaluate a gradient of a Function at the quadrature points in all cells. This evaluated gradient can then be used as input to a non-FEniCS function that calculates a material constitutive model.

Parameters:
  • e – UFL expression.

  • X – Array of points of shape (num_points, tdim) on the reference element.

  • comm – Communicator that the Expression is defined on.

  • form_compiler_options – Options used in FFCx compilation of this Expression. Run ffcx --help in the commandline to see all available options.

  • jit_options – Options controlling JIT compilation of C code.

Note

This wrapper is responsible for the FFCx compilation of the UFL Expr and attaching the correct data to the underlying C++ Expression.

X() ndarray[source]

Evaluation points on the reference cell

property argument_function_space: Optional[FunctionSpaceBase]

The argument function space if expression has argument

property code: str

C code strings

property dtype: dtype
eval(mesh: Mesh, cells: np.ndarray, values: Optional[np.ndarray] = None) np.ndarray[source]

Evaluate Expression in cells.

Parameters:
  • mesh – Mesh to evaluate Expression on.

  • cells – Cells of mesh` to evaluate Expression on.

  • values – Array to fill with evaluated values. If None, storage will be allocated. Otherwise must have shape (len(cells), num_points * value_size * num_all_argument_dofs)

Returns:

Expression evaluated at points for cells`.

property ufcx_expression

The compiled ufcx_expression object

property ufl_expression

Original UFL Expression

property value_size: int

Value size of the expression

class dolfinx.fem.Form(form, ufcx_form=None, code: Optional[str] = None)[source]

Bases: object

A finite element form

Note

Forms should normally be constructed using form() and not using this class initialiser. This class is combined with different base classes that depend on the scalar type used in the Form.

Parameters:
  • form – Compiled form object.

  • ufcx_form – UFCx form

  • code – Form C++ code

property code: Optional[str]

C code strings

property dtype: dtype

Scalar type of this form

property function_spaces: List[FunctionSpaceBase]

Function spaces on which this form is defined

property integral_types

Integral types in the form

property mesh: Mesh

Mesh on which this form is defined

property rank: int
property ufcx_form

The compiled ufcx_form object

class dolfinx.fem.Function(*args, **kw)[source]

Bases: Coefficient

A finite element function that is represented by a function space (domain, element and dofmap) and a vector holding the degrees-of-freedom

Initialize a finite element Function.

Parameters:
  • V – The function space that the Function is defined on.

  • x – Function degree-of-freedom vector. Typically required only when reading a saved Function from file.

  • name – Function name.

  • dtype – Scalar type. Is not set, the DOLFINx default scalar type is used.

collapse() Function[source]
copy() Function[source]

Create a copy of the Function. The function space is shared and the degree-of-freedom vector is copied.

property dtype: dtype
eval(x: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], cells: Union[_SupportsArray[dtype], _NestedSequence[_SupportsArray[dtype]], bool, int, float, complex, str, bytes, _NestedSequence[Union[bool, int, float, complex, str, bytes]]], u=None) ndarray[source]

Evaluate Function at points x, where x has shape (num_points, 3), and cells has shape (num_points,) and cell[i] is the index of the cell containing point x[i]. If the cell index is negative the point is ignored.

property function_space: FunctionSpaceBase

The FunctionSpace that the Function is defined on

interpolate(u: Union[Callable, Expression, Function], cells: Optional[ndarray] = None, nmm_interpolation_data=((), (), (), ())) None[source]

Interpolate an expression

Parameters:
  • u – The function, Expression or Function to interpolate.

  • cells – The cells to interpolate over. If None then all cells are interpolated over.

property name: str

Name of the Function.

split() tuple[dolfinx.fem.function.Function, ...][source]

Extract (any) sub-functions.

A sub-function can be extracted from a discrete function that is in a mixed, vector, or tensor FunctionSpace. The sub-function resides in the subspace of the mixed space.

Returns:

First level of subspaces of the function space.

sub(i: int) Function[source]

Return a sub-function.

Parameters:

i – Index of the sub-function to extract.

Note

The sub-functions are numbered i = 0, …, N-1, where N is the number of sub-spaces.

property vector

PETSc vector holding the degrees-of-freedom.

Upon first call, this function creates a PETSc Vec object that wraps the degree-of-freedom data. The Vec object is cached and the cached Vec is returned upon subsequent calls.

Note

Prefer :func`x` where possible.

property x: Vector

Vector holding the degrees-of-freedom.

dolfinx.fem.FunctionSpace(mesh: Mesh, element: Union[ufl.FiniteElementBase, ElementMetaData, Tuple[str, int, Tuple, bool]], form_compiler_options: Optional[dict[str, Any]] = None, jit_options: Optional[dict[str, Any]] = None) FunctionSpaceBase[source]

Create a finite element function space.

Deprecated since version 0.7: Use functionspace() (no caps) instead.

Parameters:
  • mesh – Mesh that space is defined on

  • element – Finite element description

  • form_compiler_options – Options passed to the form compiler

  • jit_options – Options controlling just-in-time compilation

Returns:

A function space.

class dolfinx.fem.FunctionSpaceBase(mesh: Mesh, element: ufl.FiniteElementBase, cppV: Union[_cpp.fem.FunctionSpace_float32, _cpp.fem.FunctionSpace_float64])[source]

Bases: FunctionSpace

A space on which Functions (fields) can be defined.

Create a finite element function space.

Note

This initialiser is for internal use and not normally called in user code. Use FunctionSpace() to create a function space.

Parameters:
  • mesh – Mesh that space is defined on

  • element – UFL finite element

  • cppV – Compiled C++ function space.

clone() FunctionSpaceBase[source]

Create a new FunctionSpace \(W\) which shares data with this FunctionSpace \(V\), but with a different unique integer ID.

This function is helpful for defining mixed problems and using blocked linear algebra. For example, a matrix block defined on the spaces \(V \times W\) where, \(V\) and \(W\) are defined on the same finite element and mesh can be identified as an off-diagonal block whereas the \(V \times V\) and \(V \times V\) matrices can be identified as diagonal blocks. This is relevant for the handling of boundary conditions.

Returns:

A new function space that shares data

collapse() tuple[dolfinx.fem.function.FunctionSpaceBase, numpy.ndarray][source]

Collapse a subspace and return a new function space and a map from new to old dofs.

Returns:

A new function space and the map from new to old degrees-of-freedom.

component()[source]

Return the component relative to the parent space.

contains(V) bool[source]

Check if a space is contained in, or is the same as (identity), this space.

Parameters:

V – The space to check to for inclusion.

Returns:

True if V is contained in, or is the same as, this space

property dofmap: DofMap

Degree-of-freedom map associated with the function space.

property element

Function space finite element.

property mesh: Mesh

Mesh on which the function space is defined.

property num_sub_spaces: int

Number of sub spaces.

sub(i: int) FunctionSpaceBase[source]

Return the i-th sub space.

Parameters:

i – The subspace index

Returns:

A subspace

tabulate_dof_coordinates() ndarray[Any, dtype[float64]][source]

Tabulate the coordinates of the degrees-of-freedom in the function space.

Returns:

Coordinates of the degrees-of-freedom.

Note

This method is only for elements with point evaluation degrees-of-freedom.

ufl_function_space() FunctionSpace[source]

UFL function space.

class dolfinx.fem.IntegralType(self: dolfinx.cpp.fem.IntegralType, value: int)

Bases: pybind11_object

Members:

cell

exterior_facet

interior_facet

vertex

cell = <IntegralType.cell: 0>
exterior_facet = <IntegralType.exterior_facet: 1>
interior_facet = <IntegralType.interior_facet: 2>
property name
property value
vertex = <IntegralType.vertex: 3>
dolfinx.fem.VectorFunctionSpace(mesh: Mesh, element: Union[ElementMetaData, Tuple[str, int]], dim: Optional[int] = None) FunctionSpaceBase[source]

Create a vector finite element (composition of scalar elements) function space.

Deprecated since version 0.7: Use FunctionSpace() with a shape argument instead.

Parameters:
  • mesh – Mesh that space is defined on

  • element – Finite element description. Must be a scalar element, e.g. Lagrange.

  • dim – Dimension of the vector, e.g. number of vector components. It defaults to the geometric dimension of the mesh.

Returns:

A blocked vector function space.

dolfinx.fem.apply_lifting(b: ndarray, a: List[Form], bcs: List[List[DirichletBC]], x0: Optional[List[ndarray]] = None, scale: float = 1.0, constants=None, coeffs=None) None[source]

Modify RHS vector b for lifting of Dirichlet boundary conditions.

It modifies b such that:

\[b \leftarrow b - \text{scale} * A_j (g_j - x0_j)\]

where j is a block (nest) index. For a non-blocked problem j = 0. The boundary conditions bcs are on the trial spaces V_j. The forms in [a] must have the same test space as L (from which b was built), but the trial space may differ. If x0 is not supplied, then it is treated as zero.

Note

Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End.

dolfinx.fem.assemble_matrix(a: Any, bcs: Optional[List[DirichletBC]] = None, diagonal: float = 1.0, constants=None, coeffs=None, block_mode: Optional[BlockMode] = None)[source]
dolfinx.fem.assemble_matrix(A: MatrixCSR, a: Form, bcs: Optional[List[DirichletBC]] = None, diagonal: float = 1.0, constants=None, coeffs=None) MatrixCSR

Assemble bilinear form into a matrix.

Parameters:
  • a – The bilinear form assemble.

  • bcs – Boundary conditions that affect the assembled matrix. Degrees-of-freedom constrained by a boundary condition will have their rows/columns zeroed and the value diagonal set on on the matrix diagonal.

  • constants – Constants that appear in the form. If not provided, any required constants will be computed.

  • coeffs

    Coefficients that appear in the form. If not provided,

    any required coefficients will be computed.

    block_mode: Block size mode for the returned space matrix. If

    None, default is used.

Returns:

Matrix representation of the bilinear form a.

Note

The returned matrix is not finalised, i.e. ghost values are not accumulated.

dolfinx.fem.assemble_scalar(M: Form, constants=None, coeffs=None)[source]

Assemble functional. The returned value is local and not accumulated across processes.

Parameters:
  • M – The functional to compute.

  • constants – Constants that appear in the form. If not provided, any required constants will be computed.

  • coeffs – Coefficients that appear in the form. If not provided, any required coefficients will be computed.

Returns:

The computed scalar on the calling rank.

Note

Passing constants and coefficients is a performance optimisation for when a form is assembled multiple times and when (some) constants and coefficients are unchanged.

To compute the functional value on the whole domain, the output of this function is typically summed across all MPI ranks.

dolfinx.fem.assemble_vector(L: Any, constants=None, coeffs=None)[source]
dolfinx.fem.assemble_vector(L: Form, constants=None, coeffs=None) Vector
dolfinx.fem.assemble_vector(b: ndarray, L: Form, constants=None, coeffs=None)
dolfinx.fem.bcs_by_block(spaces: Iterable[Optional[FunctionSpaceBase]], bcs: Iterable[DirichletBC]) List[List[DirichletBC]][source]

Arrange Dirichlet boundary conditions by the function space that they constrain.

Given a sequence of function spaces spaces and a sequence of DirichletBC objects bcs, return a list where the ith entry is the list of DirichletBC objects whose space is contained in space[i].

dolfinx.fem.create_matrix(a: Form, block_mode: Optional[BlockMode] = None) MatrixCSR[source]

Create a sparse matrix that is compatible with a given bilinear form. :param a: Bilinear form to assemble. :param block_mode: Block mode of the CSR matrix. If None, default :param is used.:

Returns:

Assembled sparse matrix.

dolfinx.fem.create_nonmatching_meshes_interpolation_data(*args, **kwargs)

Overloaded function.

  1. create_nonmatching_meshes_interpolation_data(mesh0: dolfinx.cpp.mesh.Mesh_float32, element0: dolfinx.cpp.fem.FiniteElement_float32, mesh1: dolfinx.cpp.mesh.Mesh_float32, padding: float) -> Tuple[List[int], List[int], List[float], List[int]]

  2. create_nonmatching_meshes_interpolation_data(geometry0: dolfinx.cpp.mesh.Geometry_float32, element0: dolfinx.cpp.fem.FiniteElement_float32, mesh1: dolfinx.cpp.mesh.Mesh_float32, cells: numpy.ndarray[numpy.int32], padding: float) -> Tuple[List[int], List[int], List[float], List[int]]

  3. create_nonmatching_meshes_interpolation_data(mesh0: dolfinx.cpp.mesh.Mesh_float64, element0: dolfinx.cpp.fem.FiniteElement_float64, mesh1: dolfinx.cpp.mesh.Mesh_float64, padding: float) -> Tuple[List[int], List[int], List[float], List[int]]

  4. create_nonmatching_meshes_interpolation_data(geometry0: dolfinx.cpp.mesh.Geometry_float64, element0: dolfinx.cpp.fem.FiniteElement_float64, mesh1: dolfinx.cpp.mesh.Mesh_float64, cells: numpy.ndarray[numpy.int32], padding: float) -> Tuple[List[int], List[int], List[float], List[int]]

dolfinx.fem.create_sparsity_pattern(a: Form)[source]

Create a sparsity pattern from a bilinear form.

Parameters:

a – Bilinear form to build a sparsity pattern for.

Returns:

Sparsity pattern for the form a.

Note

The pattern is not finalised, i.e. the caller is responsible for calling assemble on the sparsity pattern.

dolfinx.fem.create_vector(L: Form) Vector[source]

Create a Vector that is compatible with a given linear form

dolfinx.fem.dirichletbc(value: Union[Function, Constant, np.ndarray], dofs: numpy.typing.NDArray[np.int32], V: Optional[dolfinx.fem.FunctionSpaceBase] = None) DirichletBC[source]

Create a representation of Dirichlet boundary condition which is imposed on a linear system.

Parameters:
  • value – Lifted boundary values function. It must have a dtype

  • property.

  • dofs – Local indices of degrees of freedom in function space to which boundary condition applies. Expects array of size (number of dofs, 2) if function space of the problem, V, is passed. Otherwise assumes function space of the problem is the same of function space of boundary values function.

  • V – Function space of a problem to which boundary conditions are applied.

Returns:

A representation of the boundary condition for modifying linear systems.

dolfinx.fem.extract_function_spaces(forms: Union[Iterable[Form], Iterable[Iterable[Form]]], index: int = 0) Iterable[Union[None, function.FunctionSpaceBase]][source]

Extract common function spaces from an array of forms. If forms is a list of linear form, this function returns of list of the corresponding test functions. If forms is a 2D array of bilinear forms, for index=0 the list common test function space for each row is returned, and if index=1 the common trial function spaces for each column are returned.

dolfinx.fem.form(form: ~typing.Union[~ufl.form.Form, ~typing.Iterable[~ufl.form.Form]], dtype: ~typing.Union[~numpy.dtype[~typing.Any], None, ~typing.Type[~typing.Any], ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]], str, ~typing.Tuple[~typing.Any, int], ~typing.Tuple[~typing.Any, ~typing.Union[~typing.SupportsIndex, ~typing.Sequence[~typing.SupportsIndex]]], ~typing.List[~typing.Any], ~numpy._typing._dtype_like._DTypeDict, ~typing.Tuple[~typing.Any, ~typing.Any]] = <class 'numpy.float64'>, form_compiler_options: ~typing.Optional[dict] = None, jit_options: ~typing.Optional[dict] = None)[source]

Create a Form or an array of Forms.

Parameters:
  • form – A UFL form or list(s) of UFL forms.

  • dtype – Scalar type to use for the compiled form.

  • form_compiler_options – See ffcx_jit

  • jit_options – See ffcx_jit.

Returns:

Compiled finite element Form.

Note

This function is responsible for the compilation of a UFL form (using FFCx) and attaching coefficients and domains specific data to the underlying C++ form. It dynamically create a Form instance with an appropriate base class for the scalar type, e.g. _cpp.fem.Form_float64().

dolfinx.fem.functionspace(mesh: Mesh, element: Union[ufl.FiniteElementBase, ElementMetaData, Tuple[str, int, Tuple, bool]], form_compiler_options: Optional[dict[str, Any]] = None, jit_options: Optional[dict[str, Any]] = None) FunctionSpaceBase[source]

Create a finite element function space.

Parameters:
  • mesh – Mesh that space is defined on

  • element – Finite element description

  • form_compiler_options – Options passed to the form compiler

  • jit_options – Options controlling just-in-time compilation

Returns:

A function space.

dolfinx.fem.locate_dofs_geometrical(V: Union[FunctionSpaceBase, Iterable[FunctionSpaceBase]], marker: Callable) ndarray[source]

Locate degrees-of-freedom geometrically using a marker function.

Parameters:
  • V – Function space(s) in which to search for degree-of-freedom indices.

  • marker – A function that takes an array of points x with shape (gdim, num_points) and returns an array of booleans of length num_points, evaluating to True for entities whose degree-of-freedom should be returned.

Returns:

An array of degree-of-freedom indices (local to the process) for degrees-of-freedom whose coordinate evaluates to True for the marker function.

If V is a list of two function spaces, then a 2-D array of shape (number of dofs, 2) is returned.

Returned degree-of-freedom indices are unique and ordered by the first column.

dolfinx.fem.locate_dofs_topological(V: Union[FunctionSpaceBase, Iterable[FunctionSpaceBase]], entity_dim: int, entities: ndarray[Any, dtype[int32]], remote: bool = True) ndarray[source]

Locate degrees-of-freedom belonging to mesh entities topologically.

Parameters:
  • V – Function space(s) in which to search for degree-of-freedom indices.

  • entity_dim – Topological dimension of entities where degrees-of-freedom are located.

  • entities – Indices of mesh entities of dimension entity_dim where degrees-of-freedom are located.

  • remote – True to return also “remotely located” degree-of-freedom indices.

Returns:

An array of degree-of-freedom indices (local to the process) for degrees-of-freedom topologically belonging to mesh entities.

If V is a list of two function spaces, then a 2-D array of shape (number of dofs, 2) is returned.

Returned degree-of-freedom indices are unique and ordered by the first column.

dolfinx.fem.set_bc(b: ndarray, bcs: List[DirichletBC], x0: Optional[ndarray] = None, scale: float = 1.0) None[source]

Insert boundary condition values into vector.

Only local (owned) entries are set, hence communication after calling this function is not required unless ghost entries need to be updated to the boundary condition value.

dolfinx.fem.transpose_dofmap(arg0: numpy.ndarray[numpy.int32], arg1: int) dolfinx.cpp.graph.AdjacencyList_int32

Build the index to (cell, local index) map from a dofmap ((cell, local index) -> index).