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

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

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

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

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

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, mat_types, ...])

Assemble bilinear forms into matrix

assemble_vector(L[, constants, coeffs])

Assemble linear form into a new PETSc vector.

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

Assemble linear forms into a monolithic vector.

assemble_vector_nest(L[, constants, coeffs])

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

create_matrix(a[, mat_type])

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

create_matrix_block(a)

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

create_matrix_nest(a)

Create a PETSc matrix (MatNest) that is compaible with a rectangular array of bilinear forms.

create_vector(L)

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

create_vector_block(L)

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

create_vector_nest(L)

Create a PETSc netsted vector (VecNest) that is compaible with a list of linear forms.

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

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

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

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

Classes

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

Class for solving a linear variational problem of the form \(a(u, v) = L(v) \, \forall v \in V\) using PETSc as a linear algebra backend.

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

Nonlinear problem class for solving the non-linear problem \(F(u, v) = 0 \ \forall v \in V\) using PETSc as the linear algebra backend.

class dolfinx.fem.petsc.LinearProblem(a: ufl.Form, L: ufl.Form, bcs: List[DirichletBCMetaClass] = [], u: _Function = None, petsc_options={}, form_compiler_params={}, jit_params={})[source]

Bases: object

Class for solving a linear variational problem 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 – Parameters that is passed to the linear algebra backend PETSc. For available choices for the ‘petsc_options’ kwarg, see the PETSc documentation.

  • form_compiler_params – Parameters used in FFCx compilation of this form. Run ffcx --help at the commandline to see all available options.

  • jit_params – Parameters used in CFFI JIT compilation of C code generated by FFCx. See python/dolfinx/jit.py for all available parameters. Takes priority over all other parameter values.

Example:

problem = LinearProblem(a, L, [bc0, bc1], petsc_options={"ksp_type": "preonly", "pc_type": "lu"})
property A: petsc4py.PETSc.Mat

Matrix operator

property L: FormMetaClass

The compiled linear form

property a: FormMetaClass

The compiled bilinear form

property b: petsc4py.PETSc.Vec

Right-hand side vector

solve() dolfinx.fem.function.Function[source]

Solve the problem.

property solver: petsc4py.PETSc.KSP

Linear solver object

class dolfinx.fem.petsc.NonlinearProblem(F: ufl.form.Form, u: _Function, bcs: List[DirichletBCMetaClass] = [], J: ufl.form.Form = None, form_compiler_params={}, jit_params={})[source]

Bases: object

Nonlinear problem class for solving the non-linear problem \(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, \(dF/du(u) du = -F(u)\).

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_params – Parameters used in FFCx compilation of this form. Run ffcx --help at the commandline to see all available options.

  • jit_params – Parameters used in CFFI JIT compilation of C code generated by FFCx. See python/dolfinx/jit.py for all available parameters. Takes priority over all other parameter values.

Example:

problem = LinearProblem(F, u, [bc0, bc1])
F(x: petsc4py.PETSc.Vec, b: petsc4py.PETSc.Vec)[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: petsc4py.PETSc.Vec, A: petsc4py.PETSc.Mat)[source]

Assemble the Jacobian matrix.

Parameters

x – The vector containing the latest solution

property L: FormMetaClass

Compiled linear form (the residual form)

property a: FormMetaClass

Compiled bilinear form (the Jacobian form)

form(x: petsc4py.PETSc.Vec)[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: PETSc.Vec, a: List[FormMetaClass], bcs: List[List[DirichletBCMetaClass]], x0: Optional[List[PETSc.Vec]] = [], scale: float = 1.0, constants=None, coeffs=None) None[source]

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

dolfinx.fem.petsc.apply_lifting_nest(b: PETSc.Vec, a: List[List[FormMetaClass]], bcs: List[DirichletBCMetaClass], x0: Optional[PETSc.Vec] = None, scale: float = 1.0, constants=None, coeffs=None) PETSc.Vec[source]

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

dolfinx.fem.petsc.assemble_matrix(a: FormMetaClass, bcs: List[DirichletBCMetaClass] = [], diagonal: float = 1.0, constants=None, coeffs=None) PETSc.Mat[source]

Assemble bilinear form into a matrix. The returned matrix is not finalised, i.e. ghost values are not accumulated.

dolfinx.fem.petsc.assemble_matrix_block(a: List[List[FormMetaClass]], bcs: List[DirichletBCMetaClass] = [], diagonal: float = 1.0, constants=None, coeffs=None) PETSc.Mat[source]

Assemble bilinear forms into matrix

dolfinx.fem.petsc.assemble_matrix_nest(a: List[List[FormMetaClass]], bcs: List[DirichletBCMetaClass] = [], mat_types=[], diagonal: float = 1.0, constants=None, coeffs=None) PETSc.Mat[source]

Assemble bilinear forms into matrix

dolfinx.fem.petsc.assemble_vector(L: FormMetaClass, constants=None, coeffs=None) PETSc.Vec[source]

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[FormMetaClass], a: List[List[FormMetaClass]], bcs: List[DirichletBCMetaClass] = [], x0: Optional[PETSc.Vec] = None, scale: float = 1.0, constants_L=None, coeffs_L=None, constants_a=None, coeffs_a=None) PETSc.Vec[source]

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: List[FormMetaClass], constants=None, coeffs=None) 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.

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

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

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[FormMetaClass]]) PETSc.Mat[source]

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

Parameters

a – 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[FormMetaClass]]) PETSc.Mat[source]

Create a PETSc matrix (MatNest) that is compaible with a rectangular array of bilinear forms.

Parameters

a – A rectangular array of bilinear forms.

Returns

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

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

Create a PETSc vector that is compaible 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[FormMetaClass]) PETSc.Vec[source]

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

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[FormMetaClass]) PETSc.Vec[source]

Create a PETSc netsted vector (VecNest) that is compaible 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.set_bc(b: PETSc.Vec, bcs: List[DirichletBCMetaClass], x0: Optional[PETSc.Vec] = None, scale: float = 1.0) None[source]

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

dolfinx.fem.petsc.set_bc_nest(b: PETSc.Vec, bcs: List[List[DirichletBCMetaClass]], x0: Optional[PETSc.Vec] = None, scale: float = 1.0) None[source]

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