dolfin.fem

Tools for assembling and manipulating finite element forms

dolfin.fem.apply_lifting(b: petsc4py.PETSc.Vec, a: List[Union[dolfin.fem.form.Form, dolfin.cpp.fem.Form]], bcs: List[List[dolfin.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.

dolfin.fem.apply_lifting_nest(b: petsc4py.PETSc.Vec, a: List[List[Union[dolfin.fem.form.Form, dolfin.cpp.fem.Form]]], bcs: List[dolfin.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.

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

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

dolfin.fem.assemble_vector(L: Union[dolfin.fem.form.Form, dolfin.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.

dolfin.fem.assemble_vector_block(L: List[Union[dolfin.fem.form.Form, dolfin.cpp.fem.Form]], a: List[List[Union[dolfin.fem.form.Form, dolfin.cpp.fem.Form]]], bcs: List[dolfin.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 zeroed and it is not finalised, i.e. ghost values are not accumulated.

dolfin.fem.assemble_vector_nest(L: Union[dolfin.fem.form.Form, dolfin.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.

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

Assemble bilinear forms into matrix

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

Assemble bilinear forms into matrix

dolfin.fem.assemble_matrix(a: Union[dolfin.fem.form.Form, dolfin.cpp.fem.Form], bcs: List[dolfin.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.

dolfin.fem.set_bc(b: petsc4py.PETSc.Vec, bcs: List[dolfin.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.

dolfin.fem.set_bc_nest(b: petsc4py.PETSc.Vec, bcs: List[List[dolfin.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.

dolfin.fem.create_coordinate_map(o)[source]

Return a compiled UFC coordinate_mapping object

class dolfin.fem.DirichletBC(V: dolfin.function.FunctionSpace, value: Union[ufl.coefficient.Coefficient, dolfin.cpp.function.Function], domain: Union[function, List[int]], method: dolfin.cpp.fem.DirichletBC.Method = Method.topological)[source]

Bases: dolfin.cpp.fem.DirichletBC

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

class dolfin.fem.DofMap(dofmap: dolfin.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.

class dolfin.fem.Form(form: ufl.form.Form, form_compiler_parameters: dict = None)[source]

Bases: ufl.form.Form

Create dolfin Form

Parameters
  • form – Pure UFL form

  • form_compiler_parameters – Parameters used in JIT FFC compilation of this form

Note

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

class dolfin.fem.FormIntegrals

Bases: pybind11_builtins.pybind11_object

Holder for integral kernels and domains

class Type(self: dolfin.cpp.fem.FormIntegrals.Type, arg0: int) → None

Bases: pybind11_builtins.pybind11_object

Members:

cell

exterior_facet

interior_facet

property name

handle) -> str

Type

(self

dolfin.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)

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

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

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

dolfin.fem.tear(V: dolfin.function.FunctionSpace) → dolfin.function.FunctionSpace[source]

For a given function space, return the corresponding discontinuous space

dolfin.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)