dolfinx.fem

Tools for assembling and manipulating finite element forms

class dolfinx.fem.DirichletBC(value: Union[ufl.coefficient.Coefficient, dolfinx.function.Function, dolfinx.cpp.function.Function], dofs: List[int], V: dolfinx.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
cell_dofs(cell_index: int)[source]
property dof_layout
property index_map
property index_map_bs
property list
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.

class dolfinx.fem.IntegralType(self: dolfinx.cpp.fem.IntegralType, arg0: 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
vertex = <IntegralType.vertex: 3>
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]
dolfinx.fem.assemble_matrix(A: petsc4py.PETSc.Mat, a: Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form], bcs: List[dolfinx.fem.dirichletbc.DirichletBC] = [], diagonal: float = 1.0) → petsc4py.PETSc.Mat

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]
dolfinx.fem.assemble_matrix_block(A: petsc4py.PETSc.Mat, 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

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]
dolfinx.fem.assemble_matrix_nest(A: petsc4py.PETSc.Mat, 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

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]
dolfinx.fem.assemble_vector(b: petsc4py.PETSc.Vec, L: Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]) → petsc4py.PETSc.Vec

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]
dolfinx.fem.assemble_vector_block(b: petsc4py.PETSc.Vec, L: List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]], a, bcs: List[dolfinx.fem.dirichletbc.DirichletBC] = [], x0: Optional[petsc4py.PETSc.Vec] = None, scale: float = 1.0) → petsc4py.PETSc.Vec

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]
dolfinx.fem.assemble_vector_nest(b: petsc4py.PETSc.Vec, L: List[Union[dolfinx.fem.form.Form, dolfinx.cpp.fem.Form]]) → petsc4py.PETSc.Vec

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_coordinate_map(comm, o)[source]

Return a compiled UFC coordinate_mapping object

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.derivative(form: ufl.form.Form, u, du, coefficient_derivatives=None) → ufl.form.Form[source]

Compute derivative of from about u (coefficient) in the direction of du (Argument)

dolfinx.fem.increase_order(V: dolfinx.function.FunctionSpace)dolfinx.function.FunctionSpace[source]

For a given function space, return the same space, but with polynomial degree increase by 1.

dolfinx.fem.locate_dofs_geometrical(V: Iterable[Union[dolfinx.cpp.function.FunctionSpace, dolfinx.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.function.FunctionSpace, dolfinx.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.

dolfinx.fem.solve(*args, **kwargs)[source]

Solve variational problem a == L or F == 0.

The following list explains the various ways in which the solve() function can be used.

1. Solving linear variational problems

A linear variational problem a(u, v) = L(v) for all v may be solved by calling solve(a == L, u, …), where a is a bilinear form, L is a linear form, u is a Function (the solution). Optional arguments may be supplied to specify boundary conditions or solver parameters. Some examples are given below:

solve(a == L, u)
solve(a == L, u, bcs=bc)
solve(a == L, u, bcs=[bc1, bc2])

solve(a == L, u, bcs=bcs,
      solver_parameters={"linear_solver": "lu"},
      form_compiler_parameters={"optimize": True})

For available choices for the ‘solver_parameters’ kwarg, look at:

info(LinearVariationalSolver.default_parameters(), True)

2. Solving nonlinear variational problems

A nonlinear variational problem F(u; v) = 0 for all v may be solved by calling solve(F == 0, u, …), where the residual F is a linear form (linear in the test function v but possibly nonlinear in the unknown u) and u is a Function (the solution). Optional arguments may be supplied to specify boundary conditions, the Jacobian form or solver parameters. If the Jacobian is not supplied, it will be computed by automatic differentiation of the residual form. Some examples are given below:

solve(F == 0, u)
solve(F == 0, u, bcs=bc)
solve(F == 0, u, bcs=[bc1, bc2])

solve(F == 0, u, bcs, J=J,
      petsc_options={"ksp_type": "preonly", "pc_type": "lu"},
      form_compiler_parameters={"optimize": True})

For available choices for the ‘solver_parameters’ kwarg, look at:

info(NonlinearVariationalSolver.default_parameters(), True)

4. Solving linear/nonlinear variational problems adaptively

Linear and nonlinear variational problems maybe solved adaptively, with automated goal-oriented error control. The automated error control algorithm is based on adaptive mesh refinement in combination with automated generation of dual-weighted residual-based error estimates and error indicators.

An adaptive solve may be invoked by giving two additional arguments to the solve call, a numerical error tolerance and a goal functional (a Form).

M = u*dx()
tol = 1.e-6

# Linear variational problem
solve(a == L, u, bcs=bc, tol=tol, M=M)

# Nonlinear problem:
solve(F == 0, u, bcs=bc, tol=tol, M=M)
dolfinx.fem.tear(V: dolfinx.function.FunctionSpace)dolfinx.function.FunctionSpace[source]

For a given function space, return the corresponding discontinuous space