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

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

Modify the right-hand side PETSc vector b to account for constraints (Dirichlet boundary conitions).

assemble_jacobian(u, jacobian, ...)

Assemble the Jacobian and preconditioner matrices.

assemble_matrix()

Assemble a bilinear form into a matrix.

assemble_residual(u, residual, jacobian, ...)

Assemble the residual into the vector b.

assemble_vector()

Assemble linear form(s) into a new PETSc vector.

assign()

Assign Function degrees-of-freedom to a vector.

create_matrix(a[, kind])

Create a PETSc matrix that is compatible with the (sequence) of bilinear form(s).

create_vector(L[, kind])

Create a PETSc vector that is compatible with a linear form(s).

discrete_curl(space0, space1)

Assemble a discrete curl operator.

discrete_gradient(space0, space1)

Assemble a discrete gradient operator.

interpolation_matrix(V0, V1)

Assemble an interpolation operator matrix for discreye interpolation between finite element spaces.

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

Set constraint (Dirchlet boundary condition) values in an vector.

Classes

LinearProblem(a, L[, bcs, u, P, kind, ...])

Class for solving a linear variational problem.

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

Class for solving nonlinear problems with SNES.

class dolfinx.fem.petsc.LinearProblem(a: Form | Iterable[Iterable[Form]], L: Form | Iterable[Form], bcs: Iterable[DirichletBC] | None = None, u: Function | Iterable[Function] | None = None, P: Form | Iterable[Iterable[Form]] | None = None, kind: str | Iterable[Iterable[str]] | None = None, petsc_options: dict | None = None, form_compiler_options: dict | None = None, jit_options: dict | None = None, entity_maps: dict[Mesh, ndarray[Any, dtype[int32]]] | 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 or a nested sequence of bilinear forms, the left-hand side of the variational problem.

  • L – Linear UFL form or a sequence of linear forms, the right-hand side of the variational problem.

  • bcs – Sequence of Dirichlet boundary conditions to apply to the variational problem and the preconditioner matrix.

  • u – Solution function. It is created if not provided.

  • P – Bilinear UFL form or a sequence of sequence of bilinear forms, used as a preconditioner.

  • kind – The PETSc matrix and vector type. See create_matrix() for options.

  • petsc_options – Options set on the underlying PETSc KSP only. For available choices for the ‘petsc_options’ kwarg, see the PETSc KSP documentation. Options on other objects (matrices, vectors) should be set explicitly by the user.

  • form_compiler_options – Options used in FFCx compilation of all forms. 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.

  • entity_maps – If any trial functions, test functions, or coefficients in the form are not defined over the same mesh as the integration domain, entity_maps must be supplied. For each key (a mesh, different to the integration domain mesh) a map should be provided relating the entities in the integration domain mesh to the entities in the key mesh e.g. for a key-value pair (msh, emap) in entity_maps, emap[i] is the entity in msh corresponding to entity i in the integration domain mesh.

Example:

problem = LinearProblem(a, L, [bc0, bc1], petsc_options={
    "ksp_type": "preonly",
    "pc_type": "lu",
    "pc_factor_mat_solver_type": "mumps"
})

problem = LinearProblem([[a00, a01], [None, a11]], [L0, L1],
                        bcs=[bc0, bc1], u=[uh0, uh1])
property A: Mat

Left-hand side matrix.

Note

The matrix has an options prefix set.

property L: Form | Iterable[Form]

The compiled linear form representing the left-hand side.

property P_mat: Mat

Preconditioner matrix.

Note

The matrix has an options prefix set.

property a: Form | Iterable[Form]

The compiled bilinear form representing the right-hand side.

property b: Vec

Right-hand side vector.

Note

The vector has an options prefix set.

property preconditioner: Form | Iterable[Form]

The compiled bilinear form representing the preconditioner.

solve() Function | Iterable[Function][source]

Solve the problem.

property solver: KSP

The PETSc KSP solver.

Note

The KSP solver has an options prefix set.

property u: Function | Iterable[Function]

Solution function.

Note

This vector does not share memory with the solution vector x.

property x: Vec

Solution vector.

Note

This vector does not share memory with the solution Function u.

Note

The vector has an options prefix set.

class dolfinx.fem.petsc.NonlinearProblem(F: Form | Sequence[Form], u: Function | Sequence[Function], bcs: Sequence[DirichletBC] | None = None, J: Form | Sequence[Sequence[Form]] | None = None, P: Form | Sequence[Sequence[Form]] | None = None, kind: str | Iterable[Iterable[str]] | None = None, form_compiler_options: dict | None = None, jit_options: dict | None = None, petsc_options: dict | None = None, entity_maps: dict[Mesh, ndarray[Any, dtype[int32]]] | None = None)[source]

Bases: object

Class for solving nonlinear problems with SNES.

Solves problems of the form \(F_i(u, v) = 0, i=0,...N\ \forall v \in V\) where \(u=(u_0,...,u_N), v=(v_0,...,v_N)\) using PETSc SNES as the non-linear solver.

By default, the underlying SNES solver uses PETSc’s default options. To use the robust combination of LU via MUMPS with a backtracking linesearch, pass:

Example:

petsc_options = {"ksp_type": "preonly",
                 "pc_type": "lu",
                 "pc_factor_mat_solver_type": "mumps",
                 "snes_linesearch_type": "bt",
}

Note

The deprecated version of this class for use with NewtonSolver has been renamed NewtonSolverNonlinearProblem.

Parameters:
  • F – UFL form(s) representing the residual \(F_i\).

  • u – Function(s) used to define the residual and Jacobian.

  • bcs – Dirichlet boundary conditions.

  • J – UFL form(s) representing the Jacobian \(J_ij = dF_i/du_j\). If not passed, derived automatically.

  • P – UFL form(s) representing the preconditioner.

  • kind – The PETSc matrix type(s) for the Jacobian and preconditioner (MatType). See dolfinx.fem.petsc.create_matrix() for more information.

  • form_compiler_options – Options used in FFCx compilation of all forms. 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.

  • petsc_options – Options that are set on the underlying PETSc SNES object only. For available choices for the ‘petsc_options’ kwarg, see the PETSc SNES documentation. Options on other objects (matrices, vectors) should be set explicitly by the user.

  • entity_maps – If any trial functions, test functions, or coefficients in the form are not defined over the same mesh as the integration domain, entity_maps must be supplied. For each key (a mesh, different to the integration domain mesh) a map should be provided relating the entities in the integration domain mesh to the entities in the key mesh e.g. for a key-value pair (msh, emap) in entity_maps, emap[i] is the entity in msh corresponding to entity i in the integration domain mesh.

property A: Mat

Jacobian matrix.

Note

The matrix has an options prefix set.

property F: Form | Iterable[Form]

The compiled residual.

property J: Form | Iterable[Iterable[Form]]

The compiled Jacobian.

property P_mat: Vec | None

Preconditioner matrix.

Note

The matrix has an options prefix set.

property b: Vec

Residual vector.

Note

The vector has an options prefix set.

property preconditioner: Form | Iterable[Iterable[Form]] | None

The compiled preconditioner.

solve() tuple[Vec, int, int][source]

Solve the problem and update the solution in the problem instance.

Note

The user is responsible for asserting convergence of the SNES solver e.g. assert converged_reason > 0. Alternatively, pass “snes_error_if_not_converged”: True and “ksp_error_if_not_converged” : True in petsc_options.

Returns:

The solution, convergence reason and number of iterations.

property solver: SNES

The SNES solver.

Note

The SNES solver has an options prefix set.

property u: Function | Iterable[Function]

Solution function.

Note

The Function does not share memory with the solution vector x.

property x: Vec

Solution vector.

Note

The vector does not share memory with the solution Function u.

Note

The vector has an options prefix set.

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, then a is a 1D sequence. If b is blocked or a nest, then a is a 2D array of forms, with the a[i] forms used to modify the block/nest vector b[i].

  • bcs

    Boundary conditions used to modify b (see dolfinx.fem.apply_lifting()). Two cases are supported:

    1. The boundary conditions bcs are a ‘sequence-of-sequences’ such that bcs[j] are the Dirichlet boundary conditionns associated with the forms in the j th colulmn of a. Helper functions exist to create a sequence-of-sequences of DirichletBC from the 2D a and a flat Sequence of DirichletBC objects bcs:

      bcs1 = fem.bcs_by_block(
       fem.extract_function_spaces(a, 1), bcs
      )
      
    2. bcs is a sequence of dolfinx.fem.DirichletBC objects. The function deduces which DiricletBC objects apply to each column of a by matching the dolfinx.fem.FunctionSpace.

  • x0 – Vector to use in modify b (see dolfinx.fem.apply_lifting()). Treated as zero if None.

  • alpha – Scalar parameter in lifting (see dolfinx.fem.apply_lifting()).

  • constants – Packed constant data appearing in the forms a. If None, the constant data will be packed by the function.

  • coeffs – Packed coefficient data appearing in the forms a. If None, 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. Use dolfinx.fem.DirichletBC.set() to set values in b.

dolfinx.fem.petsc.assemble_jacobian(u: Sequence[Function] | Function, jacobian: Form | Iterable[Iterable[Form]], preconditioner: Form | Iterable[Iterable[Form]] | None, bcs: Iterable[DirichletBC], _snes: SNES, x: Vec, J: Mat, P: Mat)[source]

Assemble the Jacobian and preconditioner matrices.

A function conforming to the interface expected by SNES.setJacobian can be created by fixing the first four arguments:

functools.partial(assemble_jacobian, u, jacobian, preconditioner,

bcs)

Parameters:
  • u – Function tied to the solution vector within the residual and jacobian.

  • jacobian – Compiled form of the Jacobian.

  • preconditioner – Compiled form of the preconditioner.

  • bcs – List of Dirichlet boundary conditions to apply to the Jacobian and preconditioner matrices.

  • _snes – The solver instance.

  • x – The vector containing the point to evaluate at.

  • J – Matrix to assemble the Jacobian into.

  • P – Matrix to assemble the preconditioner into.

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:

  1. If a is a single bilinear form, the form is assembled into PETSc matrix of type kind.

  2. If a is a \(m imes n\) rectangular array of forms the forms in a 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 form a[i][j].

    1. If kind is a PETSc.Mat.Type (other than PETSc.Mat.Type.NEST) or is None, the matrix type is kind of the default type (if kind is None).

    2. If kind is PETSc.Mat.Type.NEST or a rectangular array of PETSc matrix types, the returned matrix has type PETSc.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() and dolfinx.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_residual(u: Function | Sequence[Function], residual: Form | Iterable[Form], jacobian: Form | Iterable[Iterable[Form]], bcs: Iterable[DirichletBC], _snes: SNES, x: Vec, b: Vec)[source]

Assemble the residual into the vector b.

A function conforming to the interface expected by SNES.setResidual can be created by fixing the first four arguments:

functools.partial(assemble_residual, u, residual, jacobian, bcs)

Parameters:
  • u – Function(s) tied to the solution vector within the residual and Jacobian.

  • residual – Form of the residual. It can be a sequence of forms.

  • jacobian – Form of the Jacobian. It can be a nested sequence of forms.

  • bcs – List of Dirichlet boundary conditions to lift the residual.

  • _snes – The solver instance.

  • x – The vector containing the point to evaluate the residual at.

  • b – Vector to assemble the residual into.

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:

  1. If L is a single linear form, the form is assembled into a ghosted PETSc vector.

  2. If L is a sequence of linear forms and kind is None or is PETSc.Vec.Type.MPI, the forms are assembled into a vector b such that b = [b_0, b_1, ..., b_n, b_0g, b_1g, ..., b_ng] where b_i are the entries associated with the ‘owned’ degrees-of-freedom for L[i] and b_ig are the ‘unowned’ (ghost) entries for L[i].

    For this case, the returned vector has an attribute _blocks that holds the local offsets into b for the (i) owned and (ii) ghost entries for each L[i]. See create_vector() for a description of the offset blocks.

  3. If L is a sequence of linear forms and kind is PETSc.Vec.Type.NEST, the forms are assembled into a PETSc nested vector b (a nest of ghosted PETSc vectors) such that L[i] is assembled into into the ith nested matrix in b.

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() and dolfinx.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 sequence of linear forms.

  • constants – Constants appearing in the form. For a single form, constants.ndim==1. For multiple forms, the constants for form L[i] are constants[i].

  • coeffs – Coefficients appearing in the form. For a single form, coeffs.shape=(num_cells, n). For multiple forms, the coefficients for form L[i] are coeffs[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 of Function``s, to ``x. When u is a sequence of Function``s, degrees-of-freedom for the ``Function``s in ``u are ‘stacked’ and assigned to x. See assign() for documentation on how stacked assignment is handled.

Parameters:
  • uFunction (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:

  1. For a single bilinear form, it creates a compatible PETSc matrix of type kind.

  2. For a rectangular array of bilinear forms, if kind is PETSc.Mat.Type.NEST or kind is an array of PETSc Mat types (with the same shape as a), a matrix of type PETSc.Mat.Type.NEST is created. The matrix is compatible with the forms a.

  3. 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 forms a. If kind is None, 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 sequence 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:

  1. For a single linear form L, if kind is None or is PETSc.Vec.Type.MPI, a ghosted PETSc vector which is compatible with L is created.

  2. If L is a sequence of linear forms and kind is None or is PETSc.Vec.Type.MPI, a ghosted PETSc vector which is compatible with L is created. The created vector b is initialized such that on each MPI process b = [b_0, b_1, ..., b_n, b_0g, b_1g, ..., b_ng], where b_i are the entries associated with the ‘owned’ degrees-of-freedom for L[i] and b_ig are the ‘unowned’ (ghost) entries for L[i].

    For this case, the returned vector has an attribute _blocks that holds the local offsets into b for the (i) owned and (ii) ghost entries for each L[i]. It can be accessed by b.getAttr("_blocks"). The offsets can be used to get views into b 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]]
    
  3. If L is a sequence of linear forms and kind is PETSc.Vec.Type.NEST, a PETSc nested vector (a ‘nest’ of ghosted PETSc vectors) which is compatible with L 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), where g 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 and bcs[i] are the boundary conditions to apply to block/nest i. Otherwise bcs should be a sequence of DirichletBCs. For block/nest problems, dolfinx.fem.bcs_by_block() can be used to prepare the 2D array of DirichletBC 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.