dolfinx.fem

Tools for assembling and manipulating finite element forms

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

Bases: ufl.constant.Constant

A constant with respect to a domain.

Parameters
  • domain (DOLFINx or UFL mesh) –

  • c – Value of the constant.

property value

Returns value of the constant.

class dolfinx.fem.DirichletBC(value: Union[ufl.coefficient.Coefficient, dolfinx.fem.function.Function, dolfinx.cpp.fem.Function], dofs: List[int], V: Optional[dolfinx.fem.function.FunctionSpace] = None)[source]

Bases: dolfinx.cpp.fem.DirichletBC

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

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 (optional) – Function space of a problem to which boundary conditions are applied.

class dolfinx.fem.DofMap(dofmap: dolfinx.cpp.fem.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 ufc_dofmap on a specific mesh.

property bs

Returns the block size of the dofmap

cell_dofs(cell_index: int)[source]

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

Parameters

cell – The cell index

property dof_layout

Returns the layout of dofs on an element

property index_map

Returns the index map that described the parallel distribution of the dofmap

property index_map_bs

Returns the block size of the index map

property list

Returns the adjacency list with dof indices for each cell

class dolfinx.fem.Expression(ufl_expression: ufl.core.expr.Expr, x: numpy.ndarray, form_compiler_parameters: dict = {}, jit_parameters: dict = {})[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
  • ufl_expression – Pure UFL expression

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

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

  • jit_parameters – Parameters 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.

property code

Return C code strings

eval(cells: numpy.ndarray, u: Optional[numpy.ndarray] = None) numpy.ndarray[source]

Evaluate Expression in cells.

Parameters
  • cells – local indices of cells to evaluate expression.

  • u (optional) – array of shape (num_cells, num_points*value_size) to store result of expression evaluation.

Returns

u – The i-th row of u contains the expression evaluated on cells[i].

Return type

np.ndarray

Note

This function allocates u of the appropriate size if u is not passed.

property num_points

Return the number of evaluation points on the reference cell.

property ufc_expression

Return the compiled ufc_expression object

property ufl_expression

Return the original UFL Expression

property value_size

Return the value size of the expression

property x

Return the evaluation points on the reference cell

class dolfinx.fem.Form(form: ufl.form.Form, form_compiler_parameters: dict = {}, jit_parameters: dict = {})[source]

Bases: object

Create DOLFINx Form

Parameters
  • form – Pure UFL form

  • form_compiler_parameters – See ffcx_jit

  • jit_parameters – See ffcx_jit

Note

This wrapper for UFL form is responsible for the actual FFCx compilation and attaching coefficients and domains specific data to the underlying C++ Form.

property code

Return C code strings

property ufc_form

Return the compiled ufc_form object

class dolfinx.fem.Function(V: dolfinx.fem.function.FunctionSpace, x: Optional[dolfinx.cpp.la.Vector] = None, name: Optional[str] = None)[source]

Bases: ufl.coefficient.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 finite element Function.

collapse()[source]
compute_point_values()[source]
copy()[source]

Return a copy of the Function. The FunctionSpace is shared and the degree-of-freedom vector is copied.

eval(x: numpy.ndarray, cells: numpy.ndarray, u=None) numpy.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: dolfinx.fem.function.FunctionSpace

Return the FunctionSpace

property id: int

Return object id index.

interpolate(u) None[source]

Interpolate an expression

property name: str

Name of the Function.

split()[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.

sub(i: int)[source]

Return a sub function.

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

ufl_evaluate(x, component, derivatives)[source]

Function used by ufl to evaluate the Expression

property vector

Return the vector holding Function degrees-of-freedom.

property x

Return the vector holding Function degrees-of-freedom.

class dolfinx.fem.FunctionSpace(mesh: dolfinx.cpp.mesh.Mesh, element: Union[ufl.finiteelement.finiteelementbase.FiniteElementBase, dolfinx.fem.function.ElementMetaData], cppV: Optional[dolfinx.cpp.fem.FunctionSpace] = None, form_compiler_parameters: dict = {}, jit_parameters: dict = {})[source]

Bases: ufl.functionspace.FunctionSpace

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

Create a finite element function space.

clone() dolfinx.fem.function.FunctionSpace[source]

Return 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.

collapse(collapsed_dofs: bool = False)[source]
Collapse a subspace and return a new function space and a map from

new to old dofs.

Arguments
collapsed_dofs

Return the map from new to old dofs

Returns
FunctionSpace

The new function space.

dict

The map from new to old dofs (optional)

component()[source]

Return the component relative to the parent space.

contains(V) bool[source]

Check whether a FunctionSpace is in this FunctionSpace, or is the same as this FunctionSpace.

property dofmap: dolfinx.fem.dofmap.DofMap

Return the degree-of-freedom map associated with the function space.

dolfin_element()[source]

Return the DOLFINx element.

property element
property id: int

The unique identifier

property mesh

Return the mesh on which the function space is defined.

num_sub_spaces() int[source]

Return the number of sub spaces.

sub(i: int) dolfinx.fem.function.FunctionSpace[source]

Return the i-th sub space.

tabulate_dof_coordinates()[source]
ufl_cell()[source]
ufl_function_space() ufl.functionspace.FunctionSpace[source]

Return the UFL function space

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

Bases: pybind11_builtins.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>
class dolfinx.fem.LinearProblem(a: ufl.form.Form, L: ufl.form.Form, bcs: List[dolfinx.fem.dirichletbc.DirichletBC] = [], u: Optional[dolfinx.fem.function.Function] = None, petsc_options={}, form_compiler_parameters={}, jit_parameters={})[source]

Bases: object

Class for solving a linear variational problem of the form a(u, v) = L(v) for all v using PETSc as a linear algebra backend.

Initialize solver for a linear variational problem.

Parameters
  • a – A bilinear UFL form, the left hand side of the variational problem.

  • L – A linear UFL form, the right hand side of the variational problem.

  • bcs – A list of Dirichlet boundary conditions.

  • u – The solution function. It will be created if not provided.

  • petsc_options – Parameters that is passed to the linear algebra backend PETSc. For available choices for the ‘petsc_options’ kwarg, see the PETSc-documentation <https://www.mcs.anl.gov/petsc/documentation/index.html>.

  • form_compiler_parameters – Parameters used in FFCx compilation of this form. Run ffcx –help at the commandline to see all available options. Takes priority over all other parameter values, except for scalar_type which is determined by DOLFINx.

  • jit_parameters – Parameters used in CFFI JIT compilation of C code generated by FFCx. See python/dolfinx/jit.py for all available parameters. Takes priority over all other parameter values.

  • code-block: (.) – python: problem = LinearProblem(a, L, [bc0, bc1], petsc_options={“ksp_type”: “preonly”, “pc_type”: “lu”})

property A: petsc4py.PETSc.Mat

Get the matrix operator

property L: dolfinx.fem.form.Form

Get the compiled linear form

property a: dolfinx.fem.form.Form

Get the compiled bilinear form

property b: petsc4py.PETSc.Vec

Get the RHS vector

solve() dolfinx.fem.function.Function[source]

Solve the problem.

property solver: petsc4py.PETSc.KSP

Get the linear solver

class dolfinx.fem.NonlinearProblem(F: ufl.form.Form, u: dolfinx.fem.function.Function, bcs: List[dolfinx.fem.dirichletbc.DirichletBC] = [], J: Optional[ufl.form.Form] = None, form_compiler_parameters={}, jit_parameters={})[source]

Bases: object

Nonlinear problem class for solving the non-linear problem F(u, v) = 0 for all v in V

Initialize class that sets up structures for solving the non-linear problem using Newton’s method, dF/du(u) du = -F(u)

Parameters
  • F – The PDE residual F(u, v)

  • u – The unknown

  • bcs – List of Dirichlet boundary conditions

  • J – UFL representation of the Jacobian (Optional)

  • form_compiler_parameters – Parameters used in FFCX compilation of this form. Run ffcx –help at the commandline to see all available options. Takes priority over all other parameter values, except for scalar_type which is determined by DOLFINX.

  • jit_parameters – Parameters used in CFFI JIT compilation of C code generated by FFCX. See python/dolfinx/jit.py for all available parameters. Takes priority over all other parameter values.

  • code-block: (.) – python: problem = LinearProblem(F, u, [bc0, bc1])

F(x: petsc4py.PETSc.Vec, b: petsc4py.PETSc.Vec)[source]

Assemble the residual F into the vector b. :param x: The vector containing the latest solution :param b: Vector to assemble the residual into

J(x: petsc4py.PETSc.Vec, A: petsc4py.PETSc.Mat)[source]

Assemble the Jacobian matrix. :param x: The vector containing the latest solution :param A: The matrix to assemble the Jacobian into

property L: dolfinx.fem.form.Form

Get the compiled linear form (the residual)

property a: dolfinx.fem.form.Form

Get the compiled bilinear form (the Jacobian)

form(x: petsc4py.PETSc.Vec)[source]

This function is called before the residual or Jacobian is computed. This is usually used to update ghost values. :param x: The vector containing the latest solution

dolfinx.fem.TensorFunctionSpace(mesh: dolfinx.cpp.mesh.Mesh, element: dolfinx.fem.function.ElementMetaData, shape=None, symmetry: Optional[bool] = None, restriction=None) dolfinx.fem.function.FunctionSpace[source]

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

dolfinx.fem.VectorFunctionSpace(mesh: dolfinx.cpp.mesh.Mesh, element: dolfinx.fem.function.ElementMetaData, dim=None, restriction=None) dolfinx.fem.function.FunctionSpace[source]

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

dolfinx.fem.adjoint(form: ufl.form.Form, reordered_arguments=None) ufl.form.Form[source]

Compute adjoint of a bilinear form by changing the ordering (count) of the test and trial functions.

The functions wraps ufl.adjoint, and by default UFL will create new Argument s. To specify the Argument s rather than creating new ones, pass a tuple of Argument s as reordered_arguments. See the documentation for ufl.adjoint for more details.

dolfinx.fem.apply_lifting(b: petsc4py.PETSc.Vec, a: List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]], bcs: List[List[dolfinx.fem.dirichletbc.DirichletBC]], x0: Optional[List[petsc4py.PETSc.Vec]] = [], scale: float = 1.0) None[source]

Modify RHS vector b for lifting of Dirichlet boundary conditions. It modifies b such that:

b <- b - 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.

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

dolfinx.fem.apply_lifting_nest(b: petsc4py.PETSc.Vec, a: List[List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]]], bcs: List[dolfinx.fem.dirichletbc.DirichletBC], x0: Optional[petsc4py.PETSc.Vec] = None, scale: float = 1.0) petsc4py.PETSc.Vec[source]

Modify nested vector for lifting of Dirichlet boundary conditions.

dolfinx.fem.assemble_matrix(a: Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form], bcs: List[dolfinx.fem.dirichletbc.DirichletBC] = [], diagonal: float = 1.0) petsc4py.PETSc.Mat[source]

Assemble bilinear form into a matrix. The returned matrix is not finalised, i.e. ghost values are not accumulated.

dolfinx.fem.assemble_matrix_block(a: List[List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]]], bcs: List[dolfinx.fem.dirichletbc.DirichletBC] = [], diagonal: float = 1.0) petsc4py.PETSc.Mat[source]

Assemble bilinear forms into matrix

dolfinx.fem.assemble_matrix_nest(a: List[List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]]], bcs: List[dolfinx.fem.dirichletbc.DirichletBC] = [], mat_types=[], diagonal: float = 1.0) petsc4py.PETSc.Mat[source]

Assemble bilinear forms into matrix

dolfinx.fem.assemble_scalar(M: Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]) numpy.float64[source]

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

dolfinx.fem.assemble_vector(L: Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]) petsc4py.PETSc.Vec[source]

Assemble linear form into a new PETSc vector. The returned vector is not finalised, i.e. ghost values are not accumulated on the owning processes.

dolfinx.fem.assemble_vector_block(L: List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]], a: List[List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]]], bcs: List[dolfinx.fem.dirichletbc.DirichletBC] = [], x0: Optional[petsc4py.PETSc.Vec] = None, scale: float = 1.0) petsc4py.PETSc.Vec[source]

Assemble linear forms into a monolithic vector. The vector is not finalised, i.e. ghost values are not accumulated.

dolfinx.fem.assemble_vector_nest(L: Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]) petsc4py.PETSc.Vec[source]

Assemble linear forms into a new nested PETSc (VecNest) vector. The returned vector is not finalised, i.e. ghost values are not accumulated on the owning processes.

dolfinx.fem.create_matrix(a: Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form], mat_type=None) petsc4py.PETSc.Mat[source]
dolfinx.fem.create_matrix_block(a: List[List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]]]) petsc4py.PETSc.Mat[source]
dolfinx.fem.create_matrix_nest(a: List[List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]]]) petsc4py.PETSc.Mat[source]
dolfinx.fem.create_vector(L: Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]) petsc4py.PETSc.Vec[source]
dolfinx.fem.create_vector_block(L: List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]]) petsc4py.PETSc.Vec[source]
dolfinx.fem.create_vector_nest(L: List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]]) petsc4py.PETSc.Vec[source]
dolfinx.fem.locate_dofs_geometrical(V: Iterable[Union[dolfinx.cpp.fem.FunctionSpace, dolfinx.fem.function.FunctionSpace]], marker: function)[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.

Return type

numpy.ndarray

dolfinx.fem.locate_dofs_topological(V: Iterable[Union[dolfinx.cpp.fem.FunctionSpace, dolfinx.fem.function.FunctionSpace]], entity_dim: int, entities: List[int], remote: bool = True)[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) – 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.

Return type

numpy.ndarray

dolfinx.fem.set_bc(b: petsc4py.PETSc.Vec, bcs: List[dolfinx.fem.dirichletbc.DirichletBC], x0: Optional[petsc4py.PETSc.Vec] = 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.set_bc_nest(b: petsc4py.PETSc.Vec, bcs: List[List[dolfinx.fem.dirichletbc.DirichletBC]], x0: Optional[petsc4py.PETSc.Vec] = None, scale: float = 1.0) None[source]

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