assemble

dolfin.fem.assembling.assemble(form, tensor=None, bcs=None, form_compiler_parameters=None, add_values=False, finalize_tensor=True, keep_diagonal=False, backend=None, cell_domains=None, exterior_facet_domains=None, interior_facet_domains=None)

Assemble the given form and return the corresponding tensor.

Arguments
Depending on the input form, which may be a functional, linear form, bilinear form or higher rank form, a scalar value, a vector, a matrix or a higher rank tensor is returned.

In the simplest case, no additional arguments are needed. However, additional arguments may and must in some cases be provided as outlined below.

The form can be either a UFL Form or a DOLFIN Form.

If the form defines integrals over different subdomains, MeshFunctions over the corresponding topological entities defining the subdomains can be provided.

The implementation of the returned tensor is determined by the default linear algebra backend. This can be overridden by specifying a different backend.

Each call to assemble() will create a new tensor. If the tensor argument is provided, this will be used instead. Sparsity pattern of tensor is reset iff tensor.empty() holds.

If the keep_diagonal is set to True, assembler ensures that potential zeros on a matrix diagonal are kept in sparsity pattern so every diagonal entry can be changed in a future (for example by ident() or ident_zeros()).

Specific form compiler parameters can be provided by the form_compiler_parameters argument. Form compiler parameters can also be controlled using the global parameters stored in parameters[“form_compiler”].

Examples of usage

The standard stiffness matrix A and a load vector b can be assembled as follows:

A = assemble(inner(grad(u),grad(v))*dx)
b = assemble(f*v*dx)

It is possible to explicitly prescribe the domains over which integrals wll be evaluated using the arguments cell_domains, exterior_facet_domains and interior_facet_domains. For instance, using a mesh function marking parts of the boundary:

# MeshFunction marking boundary parts
boundary_parts = MeshFunction("size_t", mesh, mesh.topology().dim()-1)

# Sample variational forms
a = inner(grad(u), grad(v))*dx + p*u*v*ds(0)
L = f*v*dx - g*v*ds(1) + p*q*v*ds(0)

A = assemble(a, exterior_facet_domains=boundary_parts)
b = assemble(L, exterior_facet_domains=boundary_parts)

To ensure that the assembled matrix has the right type, one may use the tensor argument:

A = PETScMatrix()
assemble(a, tensor=A)

The form a is now assembled into the PETScMatrix A.