assemble¶
-
dolfin.fem.assembling.
assemble
(form, tensor=None, form_compiler_parameters=None, add_values=False, finalize_tensor=True, keep_diagonal=False, backend=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 oftensor
is reset ifftensor.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 vectorb
can be assembled as follows:A = assemble(inner(grad(u),grad(v))*dx) b = assemble(f*v*dx)
To prescribe the domains over which integrals will be evaluated, create a Measure with the MeshFunction passed as subdomain_data. For instance, using a mesh function marking parts of the boundary:
# MeshFunction marking boundary parts boundary_markers = FacetFunction("size_t", mesh) # ... fill values in boundary_markers # Measures with references to cell and boundary markers ds = Measure("ds", subdomain_data=boundary_markers) # 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) b = assemble(L)
For interior facet integrals (dS), cell markers can be used to specify which cell is ‘+’ and which is ‘-‘. The ‘+’ and ‘-‘ sides are chosen such that the cell marker value in the cell at the ‘+’ side cell is larger than the cell marker value in the cell at the ‘-‘ side cell. If the values are equal or the cell markers are not provided, the sides are chosen arbitrarily.
Note that currently, cell markers must be associated with a cell type integral (dx), and if you don’t have such an integral a workaround is to add the integral of something over an empty domain such as ‘f*dx(99)’ with 99 a number not occuring among the cell markers. A better interface to handle this case will be provided later.
# MeshFunctions marking boundary and cell parts boundary_markers = FacetFunction("size_t", mesh) cell_markers = CellFunction("size_t", mesh) # ... fill values in boundary_markers # Measures with references to cell and boundary markers ds = Measure("ds", domain=mesh, subdomain_data=boundary_markers) dx = Measure("dx", domain=mesh, subdomain_data=cell_markers) # Sample variational forms a = inner(grad(u), grad(v))*dx + p*u*v*ds(0) L = v*dx(99) - g*v*ds(1) + p*q*v*ds(0) A = assemble(a) b = assemble(L)
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 PETScMatrixA
.