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(-> ~petsc4py.PETSc.Mat)

assemble_matrix_block(-> ~petsc4py.PETSc.Mat)

assemble_matrix_nest(-> ~petsc4py.PETSc.Mat)

assemble_vector(-> ~petsc4py.PETSc.Vec)

assemble_vector_block(-> ~petsc4py.PETSc.Vec)

assemble_vector_nest(-> ~petsc4py.PETSc.Vec)

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.DirichletBCMetaClass(value: Union[Function, Constant, numpy.ndarray], dofs: numpy.typing.ArrayLike, V: dolfinx.fem.FunctionSpace = None)[source]

Bases: object

Representation of Dirichlet boundary condition which is imposed on a linear system.

Notes

Dirichlet boundary conditions should normally be constructed using fem.dirichletbc() and not using this class initialiser. This class is combined with different base classes that depend on the scalar type of the boundary condition.

Parameters
  • value – Lifted boundary values function.

  • dofs – Local indices of degrees of freedom in function space to which boundary condition applies. Expects array of size (number of dofs, 2) if function space of the problem, V is passed. Otherwise assumes function space of the problem is the same of function space of boundary values function.

  • V – Function space of a problem to which boundary conditions are applied.

property function_space: FunctionSpace

The function space on which the boundary condition is defined

property g

The boundary condition value(s)

class dolfinx.fem.petsc.FormMetaClass(form, V: list[dolfinx.cpp.fem.FunctionSpace], coeffs, constants, subdomains: dict[dolfinx.cpp.mesh.MeshTags_int32, Optional[Any]], mesh: Mesh, ffi, code)[source]

Bases: object

A finite element form

Notes

Forms should normally be constructed using forms.form() and not using this class initialiser. This class is combined with different base classes that depend on the scalar type used in the Form.

Parameters
  • form – Compiled UFC form

  • V – The argument function spaces

  • coeffs – Finite element coefficients that appear in the form

  • constants – Constants appearing in the form

  • subdomains – Subdomains for integrals

  • mesh – The mesh that the form is deined on

property code: str

C code strings

property dtype: dtype

dtype of this form

property function_spaces: List[FunctionSpace]

Function spaces on which this form is defined

property integral_types

Integral types in the form

property mesh: Mesh

Mesh on which this form is defined

property ufcx_form

The compiled ufcx_form object

class dolfinx.fem.petsc.LinearProblem(a: Form, L: Form, bcs: List[DirichletBCMetaClass] = [], u: Optional[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: Mat

Matrix operator

property L: FormMetaClass

The compiled linear form

property a: FormMetaClass

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[DirichletBCMetaClass] = [], J: Optional[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: Vec, b: 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: Vec, A: 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: 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: Vec, a: List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]], bcs: List[List[DirichletBCMetaClass]], x0: List[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: Vec, a: List[List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]], bcs: List[DirichletBCMetaClass], x0: Optional[Vec] = None, scale: float = 1.0, 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[DirichletBCMetaClass] = [], diagonal: float = 1.0, constants=None, coeffs=None) Mat[source]
dolfinx.fem.petsc.assemble_matrix(a: Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128], bcs: List[DirichletBCMetaClass] = [], diagonal: float = 1.0, constants=None, coeffs=None) Mat
dolfinx.fem.petsc.assemble_matrix(A: Mat, a: Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128], bcs: List[DirichletBCMetaClass] = [], diagonal: float = 1.0, constants=None, coeffs=None) Mat
dolfinx.fem.petsc.assemble_matrix_block(a: Any, bcs: List[DirichletBCMetaClass] = [], diagonal: float = 1.0, constants=None, coeffs=None) Mat[source]
dolfinx.fem.petsc.assemble_matrix_block(a: List[List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]], bcs: List[DirichletBCMetaClass] = [], diagonal: float = 1.0, constants=None, coeffs=None) Mat
dolfinx.fem.petsc.assemble_matrix_block(A: Mat, a: List[List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]], bcs: List[DirichletBCMetaClass] = [], diagonal: float = 1.0, constants=None, coeffs=None) Mat
dolfinx.fem.petsc.assemble_matrix_nest(a: Any, bcs: List[DirichletBCMetaClass] = [], mat_types=[], diagonal: float = 1.0, constants=None, coeffs=None) Mat[source]
dolfinx.fem.petsc.assemble_matrix_nest(a: List[List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]], bcs: List[DirichletBCMetaClass] = [], mat_types=[], diagonal: float = 1.0, constants=None, coeffs=None) Mat
dolfinx.fem.petsc.assemble_matrix_nest(A: Mat, a: List[List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]], bcs: List[DirichletBCMetaClass] = [], diagonal: float = 1.0, constants=None, coeffs=None) Mat
dolfinx.fem.petsc.assemble_vector(L: Any, constants=None, coeffs=None) Vec[source]
dolfinx.fem.petsc.assemble_vector(L: Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128], constants=None, coeffs=None) Vec
dolfinx.fem.petsc.assemble_vector(b: Vec, L: Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128], constants=None, coeffs=None) Vec
dolfinx.fem.petsc.assemble_vector_block(L: Any, a: List[List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]], bcs: List[DirichletBCMetaClass] = [], x0: Optional[Vec] = None, scale: float = 1.0, constants_L=None, coeffs_L=None, constants_a=None, coeffs_a=None) Vec[source]
dolfinx.fem.petsc.assemble_vector_block(L: List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]], a: List[List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]], bcs: List[DirichletBCMetaClass] = [], x0: Optional[Vec] = None, scale: float = 1.0, constants_L=None, coeffs_L=None, constants_a=None, coeffs_a=None) Vec
dolfinx.fem.petsc.assemble_vector_block(b: Vec, L: List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]], a: List[List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]], bcs: List[DirichletBCMetaClass] = [], x0: Optional[Vec] = None, scale: float = 1.0, constants_L=None, coeffs_L=None, constants_a=None, coeffs_a=None) Vec
dolfinx.fem.petsc.assemble_vector_nest(L: Any, constants=None, coeffs=None) Vec[source]
dolfinx.fem.petsc.assemble_vector_nest(L: List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]], constants=None, coeffs=None) Vec
dolfinx.fem.petsc.assemble_vector_nest(b: Vec, L: List[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]], constants=None, coeffs=None) Vec
dolfinx.fem.petsc.create_matrix(a: Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128], mat_type=None) 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[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]]) 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[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]]) 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: Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]) 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[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]) 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[Union[FormMetaClass, Form_float32, Form_float64, Form_complex64, Form_complex128]]) 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: Vec, bcs: List[DirichletBCMetaClass], x0: Optional[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: Vec, bcs: List[List[DirichletBCMetaClass]], x0: Optional[Vec] = None, scale: float = 1.0) None[source]

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