dolfinx.fem.petsc
Assembly functions into PETSc objects for variational forms.
Functions in this module generally apply functions in dolfinx.fem
to PETSc linear algebra objects and handle any PETSc-specific
preparation.
Note
Due to subtle issues in the interaction between petsc4py memory
management and the Python garbage collector, it is recommended that
the PETSc method destroy()
is called on returned PETSc objects
once the object is no longer required. Note that destroy()
is
collective over the object’s MPI communicator.
Functions
|
Modify the right-hand side PETSc vector |
Assemble a bilinear form into a matrix. |
|
Assemble linear form(s) into a new PETSc vector. |
|
|
Assign |
|
Create a PETSc matrix that is compatible with the (sequence) of bilinear form(s). |
|
Create a PETSc vector that is compatible with a linear form(s). |
|
Assemble a discrete curl operator. |
|
Assemble a discrete gradient operator. |
|
Assemble an interpolation operator matrix for discreye interpolation between finite element spaces. |
|
Set constraint (Dirchlet boundary condition) values in an vector. |
Classes
|
Class for solving a linear variational problem. |
|
Nonlinear problem class for solving the non-linear problems. |
- class dolfinx.fem.petsc.LinearProblem(a: Form, L: Form, bcs: Iterable[DirichletBC] | None = None, u: Function | None = None, petsc_options: dict | None = None, form_compiler_options: dict | None = None, jit_options: dict | None = None)[source]
Bases:
object
Class for solving a linear variational problem.
Solves of the form \(a(u, v) = L(v) \, \forall v \in V\) using PETSc as a linear algebra backend.
Initialize solver for a linear variational problem.
- Parameters:
a – Bilinear UFL form, the left hand side of the variational problem.
L – Linear UFL form, the right hand side of the variational problem.
bcs – Sequecne of Dirichlet boundary conditions.
u – Solution function. It is be created if not provided.
petsc_options – Options that are passed to the linear algebra backend PETSc. For available choices for the ‘petsc_options’ kwarg, see the PETSc documentation.
form_compiler_options – Options used in FFCx compilation of this form. Run
ffcx --help
at the commandline to see all available options.jit_options – Options used in CFFI JIT compilation of C code generated by FFCx. See python/dolfinx/jit.py for all available options. Takes priority over all other option values.
Example:
problem = LinearProblem(a, L, [bc0, bc1], petsc_options={"ksp_type": "preonly", "pc_type": "lu", "pc_factor_mat_solver_type": "mumps"})
- property A: Mat
Matrix operator
- property b: Vec
Right-hand side vector
- property solver: KSP
Linear solver object
- class dolfinx.fem.petsc.NonlinearProblem(F: Form, u: Function, bcs: Iterable[DirichletBC] | None = None, J: Form = None, form_compiler_options: dict | None = None, jit_options: dict | None = None)[source]
Bases:
object
Nonlinear problem class for solving the non-linear problems.
Solves problems of the form \(F(u, v) = 0 \ \forall v \in V\) using PETSc as the linear algebra backend.
Initialize solver for solving a non-linear problem using Newton’s method`.
- Parameters:
F – The PDE residual F(u, v).
u – The unknown.
bcs – List of Dirichlet boundary conditions.
J – UFL representation of the Jacobian (optional)
form_compiler_options – Options used in FFCx compilation of this form. Run
ffcx --help
at the command line to see all available options.jit_options – Options used in CFFI JIT compilation of C code generated by FFCx. See
python/dolfinx/jit.py
for all available options. Takes priority over all other option values.
Example:
problem = NonlinearProblem(F, u, [bc0, bc1])
- F(x: Vec, b: Vec) None [source]
Assemble the residual F into the vector b.
- Parameters:
x – The vector containing the latest solution
b – Vector to assemble the residual into
- dolfinx.fem.petsc.apply_lifting(b: Vec, a: Iterable[Form] | Iterable[Iterable[Form]], bcs: Iterable[DirichletBC] | Iterable[Iterable[DirichletBC]] | None, x0: Iterable[Vec] | None = None, alpha: float = 1, constants: Iterable[ndarray] | Iterable[Iterable[ndarray]] | None = None, coeffs=None) None [source]
Modify the right-hand side PETSc vector
b
to account for constraints (Dirichlet boundary conitions).See
dolfinx.fem.apply_lifting()
for a mathematical descriptions of the lifting operation.- Parameters:
b – Vector to modify in-place.
a – List of bilinear forms. If
b
is not blocked or a nest, thena
is a 1D sequence. Ifb
is blocked or a nest, thena
is a 2D array of forms, with thea[i]
forms used to modify the block/nest vectorb[i]
.bcs –
Boundary conditions used to modify
b
(seedolfinx.fem.apply_lifting()
). Two cases are supported:The boundary conditions
bcs
is a ‘list-of-lists’ such thatbcs[j]
are the Dirichlet boundary conditionns associated with the forms in thej
th colulmn ofa
. Helper functions exist to create a list-of-lists of DirichletBC from the 2Da
and a flat list of DirichletBC objectsbcs
:bcs1 = fem.bcs_by_block(fem.extract_function_spaces(a, 1), bcs)
bcs
is a sequence ofdolfinx.fem.DirichletBC
objects. The function deduces which DiricletBC objects apply to each column ofa
by matching thedolfinx.fem.FunctionSpace
.
x0 – Vector to use in modify
b
(seedolfinx.fem.apply_lifting()
). Treated as zero ifNone
.alpha – Scalar parameter in lifting (see
dolfinx.fem.apply_lifting()
).constants – Packed constant data appearing in the forms
a
. IfNone
, the constant data will be packed by the function.coeffs – Packed coefficient data appearing in the forms
a
. IfNone
, the coefficient data will be packed by the function.
Note
Ghost contributions are not accumulated (not sent to owner). Caller is responsible for reverse-scatter to update the ghosts.
Note
Boundary condition values are not set in
b
by this function. Usedolfinx.fem.DirichletBC.set()
to set values inb
.
- dolfinx.fem.petsc.assemble_matrix(a: Form | Iterable[Iterable[Form]], bcs: Iterable[DirichletBC] | None = None, diag: float = 1, constants: Iterable[ndarray] | Iterable[Iterable[ndarray]] | None = None, coeffs: dict[tuple[_IntegralType, int], ndarray[Any, dtype[_ScalarType_co]]] | Iterable[dict[tuple[_IntegralType, int], ndarray[Any, dtype[_ScalarType_co]]]] | None = None, kind=None)[source]
- dolfinx.fem.petsc.assemble_matrix(A: Mat, a: Form | Iterable[Iterable[Form]], bcs: Iterable[DirichletBC] | None = None, diag: float = 1, constants: ndarray | Iterable[Iterable[ndarray]] | None = None, coeffs: dict[tuple[_IntegralType, int], ndarray[Any, dtype[_ScalarType_co]]] | Iterable[Iterable[dict[tuple[_IntegralType, int], ndarray[Any, dtype[_ScalarType_co]]]]] | None = None) Mat
Assemble a bilinear form into a matrix.
The following cases are supported:
If
a
is a single bilinear form, the form is assembled into PETSc matrix of typekind
.If
a
is a \(m imes n\) rectangular array of forms the forms ina
are assembled into a matrix such that:A = [A_00 ... A_0n] [A_10 ... A_1n] [ ... ] [A_m0 .. A_mn]
where
A_ij
is the matrix associated with the forma[i][j]
.If
kind
is aPETSc.Mat.Type
(other thanPETSc.Mat.Type.NEST
) or isNone
, the matrix type iskind
of the default type (ifkind
isNone
).If
kind
isPETSc.Mat.Type.NEST
or a rectangular array of PETSc matrix types, the returned matrix has typePETSc.Mat.Type.NEST
.
Rows/columns that are constrained by a Dirichlet boundary condition are zeroed, with the diagonal to set to
diag
.Constant and coefficient data that appear in the forms(s) can be packed outside of this function to avoid re-packing by this function. The functions
dolfinx.fem.pack_constants()
anddolfinx.fem.pack_coefficients()
can be used to ‘pre-pack’ the data.Note
The returned matrix is not ‘assembled’, i.e. ghost contributions are not accumulated.
- Parameters:
a – Bilinear form(s) to assembled into a matrix.
bc – Dirichlet boundary conditions applied to the system.
diag – Value to set on the matrix diagonal for Dirichlet boundary condition constrained degrees-of-freedom belonging to the same trial and test space.
constants – Constants appearing the in the form.
coeffs – Coefficients appearing the in the form.
- Returns:
Matrix representing the bilinear form.
- dolfinx.fem.petsc.assemble_vector(L: Union[Form, Iterable[Form]], constants: Optional[npt.NDArray, Iterable[npt.NDArray]] = None, coeffs: Optional[Union[dict[tuple[IntegralType, int], npt.NDArray], Iterable[dict[tuple[IntegralType, int], npt.NDArray]]]] = None, kind: Optional[str] = None) PETSc.Vec [source]
- dolfinx.fem.petsc.assemble_vector(b: PETSc.Vec, L: Union[Form, Iterable[Form]], constants: Optional[npt.NDArray, Iterable[npt.NDArray]] = None, coeffs: Optional[Union[dict[tuple[IntegralType, int], npt.NDArray], Iterable[dict[tuple[IntegralType, int], npt.NDArray]]]] = None) PETSc.Vec
Assemble linear form(s) into a new PETSc vector.
Three cases are supported:
If
L
is a single linear form, the form is assembled into a ghosted PETSc vector.If
L
is a sequence of linear forms andkind
isNone
or isPETSc.Vec.Type.MPI
, the forms are assembled into a vectorb
such thatb = [b_0, b_1, ..., b_n, b_0g, b_1g, ..., b_ng]
whereb_i
are the entries associated with the ‘owned’ degrees-of-freedom forL[i]
andb_ig
are the ‘unowned’ (ghost) entries forL[i]
.For this case, the returned vector has an attribute
_blocks
that holds the local offsets intob
for the (i) owned and (ii) ghost entries for eachL[i]
. Seecreate_vector()
for a description of the offset blocks.If
L
is a sequence of linear forms andkind
isPETSc.Vec.Type.NEST
, the forms are assembled into a PETSc nested vectorb
(a nest of ghosted PETSc vectors) such thatL[i]
is assembled into into the ith nested matrix inb
.
Constant and coefficient data that appear in the forms(s) can be packed outside of this function to avoid re-packing by this function. The functions
dolfinx.fem.pack_constants()
anddolfinx.fem.pack_coefficients()
can be used to ‘pre-pack’ the data.Note
The returned vector is not finalised, i.e. ghost values are not accumulated on the owning processes.
- Parameters:
L – A linear form or list of linear forms.
constants – Constants appearing in the form. For a single form,
constants.ndim==1
. For multiple forms, the constants for formL[i]
areconstants[i]
.coeffs – Coefficients appearing in the form. For a single form,
coeffs.shape=(num_cells, n)
. For multiple forms, the coefficients for formL[i]
arecoeffs[i]
.kind – PETSc vector type.
- Returns:
An assembled vector.
- dolfinx.fem.petsc.assign(u: Function | Sequence[Function], x: Vec)[source]
- dolfinx.fem.petsc.assign(x: Vec, u: Function | Sequence[Function])
Assign
Function
degrees-of-freedom to a vector.Assigns degree-of-freedom values in
u
, which is possibly a sequence ofFunction``s, to ``x
. Whenu
is a sequence ofFunction``s, degrees-of-freedom for the ``Function``s in ``u
are ‘stacked’ and assigned tox
. Seeassign()
for documentation on how stacked assignment is handled.- Parameters:
u –
Function
(s) to assign degree-of-freedom value from.x – Vector to assign degree-of-freedom values in
u
to.
- dolfinx.fem.petsc.create_matrix(a: Form | Iterable[Iterable[Form]], kind: str | Iterable[Iterable[str]] | None = None) Mat [source]
Create a PETSc matrix that is compatible with the (sequence) of bilinear form(s).
Three cases are supported:
For a single bilinear form, it creates a compatible PETSc matrix of type
kind
.For a rectangular array of bilinear forms, if
kind
isPETSc.Mat.Type.NEST
orkind
is an array of PETScMat
types (with the same shape asa
), a matrix of typePETSc.Mat.Type.NEST
is created. The matrix is compatible with the formsa
.For a rectangular array of bilinear forms, it create a single (non-nested) matrix of type
kind
that is compatible with the array of for formsa
. Ifkind
isNone
, then the matrix is the default type.In this case, the matrix is arranged:
A = [a_00 ... a_0n] [a_10 ... a_1n] [ ... ] [a_m0 .. a_mn]
- Parameters:
a – A bilinear form or a nested list of bilinear forms.
kind – The PETSc matrix type (
MatType
).
- Returns:
A PETSc matrix.
- dolfinx.fem.petsc.create_vector(L: Form | Iterable[Form], kind: str | None = None) Vec [source]
Create a PETSc vector that is compatible with a linear form(s).
Three cases are supported:
For a single linear form
L
, ifkind
isNone
or isPETSc.Vec.Type.MPI
, a ghosted PETSc vector which is compatible withL
is created.If
L
is a sequence of linear forms andkind
isNone
or isPETSc.Vec.Type.MPI
, a ghosted PETSc vector which is compatible withL
is created. The created vectorb
is initialized such that on each MPI processb = [b_0, b_1, ..., b_n, b_0g, b_1g, ..., b_ng]
, whereb_i
are the entries associated with the ‘owned’ degrees-of-freedom forL[i]
andb_ig
are the ‘unowned’ (ghost) entries forL[i]
.For this case, the returned vector has an attribute
_blocks
that holds the local offsets intob
for the (i) owned and (ii) ghost entries for eachL[i]
. It can be accessed byb.getAttr("_blocks")
. The offsets can be used to get views intob
for blocks, e.g.:>>> offsets0, offsets1, = b.getAttr("_blocks") >>> offsets0 (0, 12, 28) >>> offsets1 (28, 32, 35) >>> b0_owned = b.array[offsets0[0]:offsets0[1]] >>> b0_ghost = b.array[offsets1[0]:offsets1[1]] >>> b1_owned = b.array[offsets0[1]:offsets0[2]] >>> b1_ghost = b.array[offsets1[1]:offsets1[2]]
If
L
is a sequence of linear forms andkind
isPETSc.Vec.Type.NEST
, a PETSc nested vector (a ‘nest’ of ghosted PETSc vectors) which is compatible withL
is created.
- Parameters:
L – Linear form or a sequence of linear forms.
kind – PETSc vector type (
VecType
) to create.
- Returns:
A PETSc vector with a layout that is compatible with
L
. The vector is not initialised to zero.
- dolfinx.fem.petsc.discrete_curl(space0: FunctionSpace, space1: FunctionSpace) Mat [source]
Assemble a discrete curl operator.
- Parameters:
space0 – H1 space to interpolate the gradient from.
space1 – H(curl) space to interpolate into.
- Returns:
Discrete curl operator.
- dolfinx.fem.petsc.discrete_gradient(space0: FunctionSpace, space1: FunctionSpace) Mat [source]
Assemble a discrete gradient operator.
The discrete gradient operator interpolates the gradient of a H1 finite element function into a H(curl) space. It is assumed that the H1 space uses an identity map and the H(curl) space uses a covariant Piola map.
- Parameters:
space0 – H1 space to interpolate the gradient from.
space1 – H(curl) space to interpolate into.
- Returns:
Discrete gradient operator.
- dolfinx.fem.petsc.interpolation_matrix(V0: FunctionSpace, V1: FunctionSpace) Mat [source]
Assemble an interpolation operator matrix for discreye interpolation between finite element spaces.
Consider is the vector of degrees-of-freedom \(u_{i}\) associated with a function in \(V_{i}\). This function returns the matrix \(\Pi\) sucht that
\[u_{1} = \Pi u_{0}.\]- Parameters:
V0 – Space to interpolate from.
V1 – Space to interpolate into.
- Returns:
The interpolation matrix \(\Pi\).
- dolfinx.fem.petsc.set_bc(b: Vec, bcs: Iterable[DirichletBC] | Iterable[Iterable[DirichletBC]], x0: Vec | None = None, alpha: float = 1) None [source]
Set constraint (Dirchlet boundary condition) values in an vector.
For degrees-of-freedoms that are constrained by a Dirichlet boundary condition, this function sets that degrees-of-freedom to
alpha * (g - x0)
, whereg
is the boundary condition value.Only owned entries in
b
(owned by the MPI process) are modified by this function.- Parameters:
b – Vector to modify by setting boundary condition values.
bcs – Boundary conditions to apply. If
b
is nested or blocked,bcs
is a 2D array andbcs[i]
are the boundary conditions to apply to block/nesti
. Otherwisebcs
should be a sequence ofDirichletBC
s. For block/nest problems,dolfinx.fem.bcs_by_block()
can be used to prepare the 2D array ofDirichletBC
objects.x0 – Vector used in the value that constrained entries are set to. If
None
,x0
is treated as zero.alpha – Scalar value used in the value that constrained entries are set to.