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.

Functions

create_vector(L)

Create a PETSc vector that is compatible with a linear form.

create_vector_block(L)

Create a PETSc vector (blocked) that is compatible with a list of linear forms.

create_vector_nest(L)

Create a PETSc nested vector (VecNest) that is compatible with a list of linear forms.

create_matrix(a[, mat_type])

Create a PETSc matrix that is compatible with a bilinear form.

create_matrix_block(a)

Create a PETSc matrix that is compatible with a rectangular array of bilinear forms.

create_matrix_nest(a)

Create a PETSc matrix (MatNest) that is compatible with an array of bilinear forms.

assemble_vector()

Assemble linear form into a new PETSc vector.

assemble_vector_nest()

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

assemble_vector_block()

Assemble linear forms into a monolithic vector.

assemble_matrix()

Assemble bilinear form into a matrix.

assemble_matrix_nest()

Create a nested matrix and assemble bilinear forms into the matrix.

assemble_matrix_block()

Assemble bilinear forms into a blocked matrix.

apply_lifting(b, a, bcs[, x0, alpha, ...])

Apply the function dolfinx.fem.apply_lifting() to a PETSc Vector.

apply_lifting_nest(b, a, bcs[, x0, alpha, ...])

Apply the function dolfinx.fem.apply_lifting() to each sub-vector in a nested PETSc Vector.

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

Apply the function dolfinx.fem.set_bc() to a PETSc Vector.

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

Apply the function dolfinx.fem.set_bc() to each sub-vector of a nested PETSc Vector.

discrete_gradient(space0, space1)

Assemble a discrete gradient operator.

interpolation_matrix(space0, space1)

Assemble an interpolation operator matrix.

Classes

LinearProblem(a, L[, bcs, u, petsc_options, ...])

Class for solving a linear variational problem.

NonlinearProblem(F, u[, bcs, J, ...])

Nonlinear problem class for solving the non-linear problems.

class dolfinx.fem.petsc.LinearProblem(a: Form, L: Form, bcs: list[DirichletBC] = [], 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 – A bilinear UFL form, the left hand side of the variational problem.

  • L – A linear UFL form, the right hand side of the variational problem.

  • bcs – A list of Dirichlet boundary conditions.

  • u – The solution function. It will 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 L: Form

The compiled linear form

property a: Form

The compiled bilinear form

property b: Vec

Right-hand side vector

solve() Function[source]

Solve the problem.

property solver: KSP

Linear solver object

class dolfinx.fem.petsc.NonlinearProblem(F: Form, u: Function, bcs: list[DirichletBC] = [], 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 = LinearProblem(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

J(x: Vec, A: Mat) None[source]

Assemble the Jacobian matrix.

Parameters:

x – The vector containing the latest solution

property L: Form

Compiled linear form (the residual form)

property a: Form

Compiled bilinear form (the Jacobian form)

form(x: Vec) None[source]

This function is called before the residual or Jacobian is computed. This is usually used to update ghost values.

Parameters:

x – The vector containing the latest solution

dolfinx.fem.petsc.apply_lifting(b: Vec, a: list[Form], bcs: list[list[DirichletBC]], x0: list[Vec] = [], alpha: float = 1, constants=None, coeffs=None) None[source]

Apply the function dolfinx.fem.apply_lifting() to a PETSc Vector.

dolfinx.fem.petsc.apply_lifting_nest(b: Vec, a: list[list[Form]], bcs: list[DirichletBC], x0: Vec | None = None, alpha: float = 1, constants=None, coeffs=None) Vec[source]

Apply the function dolfinx.fem.apply_lifting() to each sub-vector in a nested PETSc Vector.

dolfinx.fem.petsc.assemble_matrix(a: Any, bcs: list[DirichletBC] = [], diagonal: float = 1.0, constants=None, coeffs=None)[source]
dolfinx.fem.petsc.assemble_matrix(A: Mat, a: Form, bcs: list[DirichletBC] = [], diagonal: float = 1.0, constants=None, coeffs=None) Mat

Assemble bilinear form into a matrix.

The returned matrix is not finalised, i.e. ghost values are not accumulated.

Note

The returned matrix is not ‘assembled’, i.e. ghost contributions have not been communicated.

Parameters:
  • a – Bilinear form to assembled into a matrix.

  • bc – Dirichlet boundary conditions applied to the system.

  • diagonal – 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_matrix_block(a: list[list[Form]], bcs: list[DirichletBC] = [], diagonal: float = 1.0, constants=None, coeffs=None) Mat[source]
dolfinx.fem.petsc.assemble_matrix_block(A: Mat, a: list[list[Form]], bcs: list[DirichletBC] = [], diagonal: float = 1.0, constants=None, coeffs=None) Mat

Assemble bilinear forms into a blocked matrix.

dolfinx.fem.petsc.assemble_matrix_nest(a: list[list[Form]], bcs: list[DirichletBC] = [], mat_types=[], diagonal: float = 1.0, constants=None, coeffs=None) Mat[source]
dolfinx.fem.petsc.assemble_matrix_nest(A: Mat, a: list[list[Form]], bcs: list[DirichletBC] = [], diagonal: float = 1.0, constants=None, coeffs=None) Mat

Create a nested matrix and assemble bilinear forms into the matrix.

Parameters:
  • a – Rectangular (list-of-lists) array for bilinear forms.

  • bcs – Dirichlet boundary conditions.

  • mat_types – PETSc matrix type for each matrix block.

  • diagonal – 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:

PETSc matrix (MatNest) representing the block of bilinear forms.

dolfinx.fem.petsc.assemble_vector(L: Any, constants=None, coeffs=None) Vec[source]
dolfinx.fem.petsc.assemble_vector(b: Vec, L: Form, constants=None, coeffs=None) Vec

Assemble linear form into a new PETSc vector.

Note

The returned vector is not finalised, i.e. ghost values are not accumulated on the owning processes.

Parameters:

L – A linear form.

Returns:

An assembled vector.

dolfinx.fem.petsc.assemble_vector_block(L: list[Form], a: list[list[Form]], bcs: list[DirichletBC] = [], x0: Vec | None = None, alpha: float = 1, constants_L=None, coeffs_L=None, constants_a=None, coeffs_a=None) Vec[source]
dolfinx.fem.petsc.assemble_vector_block(b: Vec, L: list[Form], a: list[list[Form]], bcs: list[DirichletBC] = [], x0: Vec | None = None, alpha: float = 1, constants_L=None, coeffs_L=None, constants_a=None, coeffs_a=None) Vec

Assemble linear forms into a monolithic vector.

The vector is not finalised, i.e. ghost values are not accumulated.

dolfinx.fem.petsc.assemble_vector_nest(L: Any, constants=None, coeffs=None) Vec[source]
dolfinx.fem.petsc.assemble_vector_nest(b: Vec, L: list[Form], constants=None, coeffs=None) Vec

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.

dolfinx.fem.petsc.create_matrix(a: Form, mat_type=None) Mat[source]

Create a PETSc matrix that is compatible with a bilinear form.

Note

Due to subtle issues in the interaction between petsc4py memory management and the Python garbage collector, it is recommended that the method PETSc.Mat.destroy() is called on the returned object once the object is no longer required. Note that PETSc.Mat.destroy() is collective over the object’s MPI communicator.

Parameters:
  • a – A bilinear form.

  • mat_type – The PETSc matrix type (MatType).

Returns:

A PETSc matrix with a layout that is compatible with a.

dolfinx.fem.petsc.create_matrix_block(a: list[list[Form]]) Mat[source]

Create a PETSc matrix that is compatible with a rectangular array of bilinear forms.

Note

Due to subtle issues in the interaction between petsc4py memory management and the Python garbage collector, it is recommended that the method PETSc.Mat.destroy() is called on the returned object once the object is no longer required. Note that PETSc.Mat.destroy() is collective over the object’s MPI communicator.

Parameters:

a – Rectangular array of bilinear forms.

Returns:

A PETSc matrix with a blocked layout that is compatible with a.

dolfinx.fem.petsc.create_matrix_nest(a: list[list[Form]]) Mat[source]

Create a PETSc matrix (MatNest) that is compatible with an array of bilinear forms.

Note

Due to subtle issues in the interaction between petsc4py memory management and the Python garbage collector, it is recommended that the method PETSc.Mat.destroy() is called on the returned object once the object is no longer required. Note that PETSc.Mat.destroy() is collective over the object’s MPI communicator.

Parameters:

a – Rectangular array of bilinear forms.

Returns:

A PETSc matrix (MatNest) that is compatible with a.

dolfinx.fem.petsc.create_vector(L: Form) Vec[source]

Create a PETSc vector that is compatible with a linear form.

Parameters:

L – A linear form.

Returns:

A PETSc vector with a layout that is compatible with L.

dolfinx.fem.petsc.create_vector_block(L: list[Form]) Vec[source]

Create a PETSc vector (blocked) that is compatible with a list of linear forms.

Note

Due to subtle issues in the interaction between petsc4py memory management and the Python garbage collector, it is recommended that the method PETSc.Vec.destroy() is called on the returned object once the object is no longer required. Note that PETSc.Vec.destroy() is collective over the object’s MPI communicator.

Parameters:

L – List of linear forms.

Returns:

A PETSc vector with a layout that is compatible with L.

dolfinx.fem.petsc.create_vector_nest(L: list[Form]) Vec[source]

Create a PETSc nested vector (VecNest) that is compatible with a list of linear forms.

Parameters:

L – List of linear forms.

Returns:

A PETSc nested vector (VecNest) with a layout that is compatible with L.

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.

Note

Due to subtle issues in the interaction between petsc4py memory management and the Python garbage collector, it is recommended that the method PETSc.Mat.destroy() is called on the returned object once the object is no longer required. Note that PETSc.Mat.destroy() is collective over the object’s MPI communicator.

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(space0: FunctionSpace, space1: FunctionSpace) Mat[source]

Assemble an interpolation operator matrix.

Note

Due to subtle issues in the interaction between petsc4py memory management and the Python garbage collector, it is recommended that the method PETSc.Mat.destroy() is called on the returned object once the object is no longer required. Note that PETSc.Mat.destroy() is collective over the object’s MPI communicator.

Parameters:
  • space0 – Space to interpolate from.

  • space1 – Space to interpolate into.

Returns:

Interpolation matrix.

dolfinx.fem.petsc.set_bc(b: Vec, bcs: list[DirichletBC], x0: Vec | None = None, alpha: float = 1) None[source]

Apply the function dolfinx.fem.set_bc() to a PETSc Vector.

dolfinx.fem.petsc.set_bc_nest(b: Vec, bcs: list[list[DirichletBC]], x0: Vec | None = None, alpha: float = 1) None[source]

Apply the function dolfinx.fem.set_bc() to each sub-vector of a nested PETSc Vector.