dolfin.fem.assemble

Assembly functions for variational forms.

Functions

apply_lifting(b, a, bcs[, x0, scale])

Modify RHS vector b for lifting of Dirichlet boundary conditions.

apply_lifting_nest(b, a, bcs[, x0, scale])

Modify nested vector for lifting of Dirichlet boundary conditions.

assemble_matrix(a[, bcs, diagonal])

Assemble bilinear form into a matrix.

assemble_matrix_block(a, bcs[, diagonal])

Assemble bilinear forms into matrix

assemble_matrix_nest(a, bcs[, diagonal])

Assemble bilinear forms into matrix

assemble_scalar(M)

Assemble functional.

assemble_vector(L)

Assemble linear form into a new PETSc vector.

assemble_vector_block(L, a, bcs[, x0, scale])

Assemble linear forms into a monolithic vector.

assemble_vector_nest(L)

Assemble linear forms into a new nested PETSc (VecNest) vector.

create_matrix(a)

create_matrix_block(a)

create_matrix_nest(a)

create_vector(L)

create_vector_block(L)

create_vector_nest(L)

set_bc(b, bcs[, x0, scale])

Insert boundary condition values into vector.

set_bc_nest(b, bcs[, x0, scale])

Insert boundary condition values into nested vector.

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