dolfinx.fem.petsc
High-level solver classes and functions for assembling PETSc objects.
Functions in this module generally apply functions in dolfinx.fem
to PETSc linear algebra objects and handle any PETSc-specific
preparation.
Note
The following does not apply to the high-level classes
dolfinx.fem.petsc.LinearProblem
dolfinx.fem.petsc.NonlinearProblem.
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 the Jacobian and preconditioner matrices at  | 
| Assemble a bilinear form into a matrix. | |
| 
 | Assemble the residual at  | 
| 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) or functionspace(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
| 
 | High-level class for solving a linear variational problem using a PETSc KSP. | 
| 
 | (Deprecated) Nonlinear problem class for solving nonlinear problems using  | 
| 
 | High-level class for solving nonlinear variational problems with PETSc SNES. | 
| Utility attributes for working with CFFI (ABI mode) and Numba. | |
| Utility attributes for working with ctypes and PETSc. | |
| Utility attributes for working with Numba and PETSc. | 
- class dolfinx.fem.petsc.LinearProblem(a: Form | Sequence[Sequence[Form]], L: Form | Sequence[Form], *, petsc_options_prefix: str, bcs: Sequence[DirichletBC] | None = None, u: Function | Sequence[Function] | None = None, P: Form | Sequence[Sequence[Form]] | None = None, kind: str | Sequence[Sequence[str]] | None = None, petsc_options: dict | None = None, form_compiler_options: dict | None = None, jit_options: dict | None = None, entity_maps: Sequence[EntityMap] | None = None)[source]
- Bases: - object- High-level class for solving a linear variational problem using a PETSc KSP. - Solves problems of the form \(a_{ij}(u, v) = f_i(v), i,j=0,\ldots,N\ \forall v \in V\) where \(u=(u_0,\ldots,u_N), v=(v_0,\ldots,v_N)\) using PETSc KSP as the linear solver. - Note - This high-level class automatically handles PETSc memory management. The user does not need to manually call - .destroy()on returned PETSc objects.- Initialize solver for a linear variational problem. - By default, the underlying KSP solver uses PETSc’s default options, usually GMRES + ILU preconditioning. To use the robust combination of LU via MUMPS - Example: - problem = LinearProblem(a, L, bcs=[bc0, bc1], petsc_options_prefix="basic_linear_problem", petsc_options= { "ksp_type": "preonly", "pc_type": "lu", "pc_factor_mat_solver_type": "mumps" }) - This class also supports nested block-structured problems. - Example: - problem = LinearProblem([[a00, a01], [None, a11]], [L0, L1], bcs=[bc0, bc1], u=[uh0, uh1], kind="nest", petsc_options_prefix="nest_linear_problem") - Every PETSc object created will have a unique options prefix set. We recommend discovering these prefixes dynamically via the petsc4py API rather than hard-coding each prefix value into the programme. - Example: - ksp_options_prefix = problem.solver.getOptionsPrefix() A_options_prefix = problem.A.getOptionsPrefix() - 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 kind. Common choices are - mpiand- nest. See- dolfinx.fem.petsc.create_matrix()and- dolfinx.fem.petsc.create_vector()for more information.
- petsc_options_prefix – Mandatory named argument. Options prefix used as root prefix on all internally created PETSc objects. Typically ends with - _. Must be the same on all ranks, and is usually unique within the programme.
- petsc_options – Options set on the underlying PETSc KSP only. The options must be the same on all ranks. For available choices for the - petsc_optionskwarg, 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 --helpat 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.pyfor 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, a corresponding :class: EntityMap<dolfinx.mesh.EntityMap> must be provided. 
 
 - property A: Mat
- Left-hand side matrix. 
 - property P_mat: Mat
- Preconditioner matrix. 
 - property b: Vec
- Right-hand side vector. 
 - property preconditioner: Form | Sequence[Form]
- The compiled bilinear form representing the preconditioner. 
 - solve() Function | Sequence[Function][source]
- Solve the problem. - This method updates the solution - ufunction(s) stored in the problem instance.- Note - The user is responsible for asserting convergence of the KSP solver e.g. - problem.solver.getConvergedReason() > 0. Alternatively, pass- "ksp_error_if_not_converged" : Truein- petsc_optionsto raise a- PETScErroron failure.- Returns:
- The solution function(s). 
 
 - property solver: KSP
- The PETSc KSP solver. 
 - property u: Function | Sequence[Function]
- Solution function(s). - Note - The function(s) do not share memory with the solution vector - x.
 - property x: Vec
- Solution vector. - Note - The vector does not share memory with the solution function(s) - u.
 
- class dolfinx.fem.petsc.NewtonSolverNonlinearProblem(F: Form, u: Function, bcs: Sequence[DirichletBC] | None = None, J: Form = None, form_compiler_options: dict | None = None, jit_options: dict | None = None)[source]
- Bases: - object- (Deprecated) Nonlinear problem class for solving nonlinear problems using - dolfinx.nls.petsc.NewtonSolver.- Solves problems of the form \(F(u, v) = 0 \ \forall v \in V\) using PETSc as the linear algebra backend. - Note - This class is deprecated in favour of - dolfinx.fem.petsc.NonlinearProblem, a high level interface to SNES.- Note - This class was previously called - dolfinx.fem.petsc.NonlinearProblem.- 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 --helpat 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.pyfor 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 
 
 
 
- class dolfinx.fem.petsc.NonlinearProblem(F: Form | Sequence[Form], u: Function | Sequence[Function], *, petsc_options_prefix: str, bcs: Sequence[DirichletBC] | None = None, J: Form | Sequence[Sequence[Form]] | None = None, P: Form | Sequence[Sequence[Form]] | None = None, kind: str | Sequence[Sequence[str]] | None = None, petsc_options: dict | None = None, form_compiler_options: dict | None = None, jit_options: dict | None = None, entity_maps: Sequence[EntityMap] | None = None)[source]
- Bases: - object- High-level class for solving nonlinear variational problems with PETSc SNES. - Solves problems of the form \(F_i(u, v) = 0, i=0,\ldots,N\ \forall v \in V\) where \(u=(u_0,\ldots,u_N), v=(v_0,\ldots,v_N)\) using PETSc SNES as the non-linear solver. - Note - The deprecated version of this class for use with - dolfinx.nls.petsc.NewtonSolverhas been renamed- dolfinx.fem.petsc.NewtonSolverNonlinearProblem.- Note - This high-level class automatically handles PETSc memory management. The user does not need to manually call - .destroy()on returned PETSc objects.- Initialize solver for a nonlinear variational problem. - 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", } - Every PETSc object will have a unique options prefix set. We recommend discovering these prefixes dynamically via the petsc4py API rather than hard-coding each prefix value into the programme. - Example: - snes_options_prefix = problem.solver.getOptionsPrefix() jacobian_options_prefix = problem.A.getOptionsPrefix() - 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 and vector kind. Common choices are - mpiand- nest. See- dolfinx.fem.petsc.create_matrix()and- dolfinx.fem.petsc.create_vector()for more information.
- petsc_options_prefix – Mandatory named argument. Options prefix used as root prefix on all internally created PETSc objects. Typically ends with _. Must be the same on all ranks, and is usually unique within the programme. 
- petsc_options – Options set on the underlying PETSc SNES only. The options must be the same on all ranks. For available choices for - petsc_options, see the PETSc SNES 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 --helpat 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.pyfor 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, a corresponding :class: EntityMap<dolfinx.mesh.EntityMap> must be provided. 
 
 - property A: Mat
- Jacobian matrix. 
 - property P_mat: Mat | None
- Preconditioner matrix. 
 - property b: Vec
- Residual vector. 
 - solve() Function | Sequence[Function][source]
- Solve the problem. - This method updates the solution - ufunction(s) stored in the problem instance.- Note - The user is responsible for asserting convergence of the SNES solver e.g. - assert problem.solver.getConvergedReason() > 0. Alternatively, pass- "snes_error_if_not_converged": Trueand- "ksp_error_if_not_converged" : Truein- petsc_optionsto raise a- PETScErroron failure.- Returns:
- The solution function(s). 
 
 - property solver: SNES
- The SNES solver. 
 - property u: Function | Sequence[Function]
- Solution function(s). - Note - The function(s) do not share memory with the solution vector - x.
 - property x: Vec
- Solution vector. - Note - The vector does not share memory with the solution function(s) - u.
 
- dolfinx.fem.petsc.apply_lifting(b: Vec, a: Sequence[Form] | Sequence[Sequence[Form]], bcs: Sequence[DirichletBC] | Sequence[Sequence[DirichletBC]] | None, x0: Sequence[Vec] | None = None, alpha: float = 1, constants: Sequence[ndarray[tuple[int, ...], dtype[_ScalarType_co]]] | Sequence[Sequence[ndarray[tuple[int, ...], dtype[_ScalarType_co]]]] | None = None, coeffs: dict[tuple[_IntegralType, int], ndarray[tuple[int, ...], dtype[_ScalarType_co]]] | Sequence[Sequence[dict[tuple[_IntegralType, int], ndarray[tuple[int, ...], dtype[_ScalarType_co]]]]] | None = None) None[source]
- Modify the right-hand side PETSc vector - bto 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 - bis not blocked or a nest, then- ais a 1D sequence. If- bis blocked or a nest, then- ais a 2D array of forms, with the- a[i]forms used to modify the block/nest vector- b[i].
- bcs – - Boundary conditions to apply, which form a 2D array. If - bis nested or blocked then- bcs[i]are the boundary conditions to apply to block/nest- i. The function- dolfinx.fem.bcs_by_block()can be used to prepare the 2D array of- DirichletBCobjects from the 2D sequence- a:- bcs1 = fem.bcs_by_block( fem.extract_function_spaces(a, 1), bcs ) - If - bis not blocked or nest, then- len(bcs)must be equal to 1. The function- dolfinx.fem.bcs_by_block()can be used to prepare the 2D array of- DirichletBCfrom the 1D sequence- a:- bcs1 = fem.bcs_by_block( fem.extract_function_spaces([a], 1), bcs ) 
- 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 - bby this function. Use- dolfinx.fem.DirichletBC.set()to set values in- b.
- dolfinx.fem.petsc.assemble_jacobian(u: Sequence[Function] | Function, jacobian: Form | Sequence[Sequence[Form]], preconditioner: Form | Sequence[Sequence[Form]] | None, bcs: Sequence[DirichletBC], _snes: SNES, x: Vec, J: Mat, P_mat: Mat)[source]
- Assemble the Jacobian and preconditioner matrices at - xinto- Jand- P_mat.- A function conforming to the interface expected by - SNES.setJacobiancan be created by fixing the first four arguments e.g.:- Example: - snes = PETSc.SNES().create(mesh.comm) assemble_jacobian = functools.partial( dolfinx.fem.petsc.assemble_jacobian, u, jacobian, preconditioner, bcs) snes.setJacobian(assemble_jacobian, A, P_mat) - 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_mat – Matrix to assemble the preconditioner into. 
 
 
- dolfinx.fem.petsc.assemble_matrix(a: Form | Sequence[Sequence[Form]], bcs: Sequence[DirichletBC] | None = None, diag: float = 1, constants: Sequence[ndarray[tuple[int, ...], dtype[_ScalarType_co]]] | Sequence[Sequence[ndarray[tuple[int, ...], dtype[_ScalarType_co]]]] | None = None, coeffs: dict[tuple[_IntegralType, int], ndarray[tuple[int, ...], dtype[_ScalarType_co]]] | Sequence[dict[tuple[_IntegralType, int], ndarray[tuple[int, ...], dtype[_ScalarType_co]]]] | None = None, kind=None)[source]
- dolfinx.fem.petsc.assemble_matrix(A: Mat, a: Form | Sequence[Sequence[Form]], bcs: Sequence[DirichletBC] | None = None, diag: float = 1, constants: ndarray[tuple[int, ...], dtype[_ScalarType_co]] | Sequence[Sequence[ndarray[tuple[int, ...], dtype[_ScalarType_co]]]] | None = None, coeffs: dict[tuple[_IntegralType, int], ndarray[tuple[int, ...], dtype[_ScalarType_co]]] | Sequence[Sequence[dict[tuple[_IntegralType, int], ndarray[tuple[int, ...], dtype[_ScalarType_co]]]]] | None = None) Mat
- Assemble a bilinear form into a matrix. - The following cases are supported: - If - ais a single bilinear form, the form is assembled into PETSc matrix of type- kind.
- If - ais a \(m \times n\) rectangular array of forms the forms in- aare assembled into a matrix such that:- A = [A_00 ... A_0n] [A_10 ... A_1n] [ ... ] [A_m0 .. A_mn] - where - A_ijis the matrix associated with the form- a[i][j].- If - kindis a- PETSc.Mat.Type(other than- PETSc.Mat.Type.NEST) or is- None, the matrix type is- kindof the default type (if- kindis- None).
- If - kindis- PETSc.Mat.Type.NESTor 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 | Sequence[Form], jacobian: Form | Sequence[Sequence[Form]], bcs: Sequence[DirichletBC], _snes: SNES, x: Vec, b: Vec)[source]
- Assemble the residual at - xinto the vector- b.- A function conforming to the interface expected by - SNES.setFunctioncan be created by fixing the first four arguments, e.g.:- Example: - snes = PETSc.SNES().create(mesh.comm) assemble_residual = functools.partial( dolfinx.fem.petsc.assemble_residual, u, residual, jacobian, bcs) snes.setFunction(assemble_residual, x, b) - 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: Form | Sequence[Form], constants: ndarray[tuple[int, ...], dtype[_ScalarType_co]] | Sequence[ndarray[tuple[int, ...], dtype[_ScalarType_co]]] | None = None, coeffs: dict[tuple[_IntegralType, int], ndarray[tuple[int, ...], dtype[_ScalarType_co]]] | Sequence[dict[tuple[_IntegralType, int], ndarray[tuple[int, ...], dtype[_ScalarType_co]]]] | None = None, kind: str | None = None) Vec[source]
- dolfinx.fem.petsc.assemble_vector(b: Vec, L: Form | Sequence[Form], constants: ndarray[tuple[int, ...], dtype[_ScalarType_co]] | Sequence[ndarray[tuple[int, ...], dtype[_ScalarType_co]]] | None = None, coeffs: dict[tuple[_IntegralType, int], ndarray[tuple[int, ...], dtype[_ScalarType_co]]] | Sequence[dict[tuple[_IntegralType, int], ndarray[tuple[int, ...], dtype[_ScalarType_co]]]] | None = None) Vec
- Assemble linear form(s) into a new PETSc vector. - Three cases are supported: - If - Lis a single linear form, the form is assembled into a ghosted PETSc vector.
- If - Lis a sequence of linear forms and- kindis- Noneor is- PETSc.Vec.Type.MPI, the forms are assembled into a vector- bsuch that- b = [b_0, b_1, ..., b_n, b_0g, b_1g, ..., b_ng]where- b_iare the entries associated with the ‘owned’ degrees-of-freedom for- L[i]and- b_igare the ‘unowned’ (ghost) entries for- L[i].- For this case, the returned vector has an attribute - _blocksthat holds the local offsets into- bfor the (i) owned and (ii) ghost entries for each- L[i]. See- create_vector()for a description of the offset blocks.
- If - Lis a sequence of linear forms and- kindis- 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 - Functiondegrees-of-freedom to a vector.- Assigns degree-of-freedom values in - u, which is possibly a sequence of- Function``s, to ``x. When- uis a sequence of- Function``s, degrees-of-freedom for the ``Function``s in ``uare ‘stacked’ and assigned to- x. See- assign()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 - uto.
 
 
- class dolfinx.fem.petsc.cffi_utils[source]
- Bases: - object- Utility attributes for working with CFFI (ABI mode) and Numba. - Registers Numba’s complex types with CFFI. - If PETSc is available, CFFI convenience functions for calling PETSc C functions are also created. These are typically called from within Numba functions. - Examples - A typical use of these utility functions is: - import numpy as np import numpy.typing as npt def set_vals(A: int, m: int, rows: npt.NDArray[PETSc.IntType], n: int, cols: npt.NDArray[PETSc.IntType], data: npt.NDArray[PETSc.ScalarType], mode: int): MatSetValuesLocal(A, m, ffi.from_buffer(rows), n, ffi.from_buffer(cols), ffi.from_buffer(rows(data), mode) - MatSetValuesBlockedLocal = <cdata 'int(*)(void *, int32_t, int32_t *, int32_t, int32_t *, double *, int)' 0x7f740c1e0760>
- See PETSc MatSetValuesBlockedLocal documentation. 
 - MatSetValuesLocal = <cdata 'int(*)(void *, int32_t, int32_t *, int32_t, int32_t *, double *, int)' 0x7f740c1c4550>
- See PETSc MatSetValuesLocal documentation. 
 
- dolfinx.fem.petsc.create_matrix(a: Form | Sequence[Sequence[Form]], kind: str | Sequence[Sequence[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 - kindis- PETSc.Mat.Type.NESTor- kindis an array of PETSc- Mattypes (with the same shape as- a), a matrix of type- PETSc.Mat.Type.NESTis created. The matrix is compatible with the forms- a.
- For a rectangular array of bilinear forms, it create a single (non-nested) matrix of type - kindthat is compatible with the array of for forms- a. If- kindis- 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(V: FunctionSpace | Sequence[FunctionSpace | None], /, kind: str | None = None) Vec[source]
- Create a PETSc vector that is compatible with a linear form(s) or functionspace(s). Three cases are supported: - For a single space - V, if- kindis- Noneor is- PETSc.Vec.Type.MPI, a ghosted PETSc vector which is compatible with- Vis created.
- If - Vis a sequence of functionspaces and- kindis- Noneor is- PETSc.Vec.Type.MPI, a ghosted PETSc vector which is compatible with- Vis created. The created vector- bis initialized such that on each MPI process- b = [b_0, b_1, ..., b_n, b_0g, b_1g, ..., b_ng], where- b_iare the entries associated with the ‘owned’ degrees-of-freedom for- V[i]and- b_igare the ‘unowned’ (ghost) entries for- V[i].- For this case, the returned vector has an attribute - _blocksthat holds the local offsets into- bfor the (i) owned and (ii) ghost entries for each- V_i. It can be accessed by- b.getAttr("_blocks"). The offsets can be used to get views into- bfor 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/Vis a sequence of linear forms/functionspaces and- kindis- PETSc.Vec.Type.NEST, a PETSc nested vector (a ‘nest’ of ghosted PETSc vectors) which is compatible with- L/Vis created.
 - Parameters:
- V – Function space or a sequence of such. 
- kind – PETSc vector type ( - VecType) to create.
 
- Returns:
- A PETSc vector with a layout that is compatible with - V. The vector is not initialised to zero.
 
- class dolfinx.fem.petsc.ctypes_utils[source]
- Bases: - object- Utility attributes for working with ctypes and PETSc. - These attributes are convenience functions for calling PETSc C functions, typically from within Numba functions. - Examples - A typical use of these utility functions is: - import numpy as np import numpy.typing as npt def set_vals(A: int, m: int, rows: npt.NDArray[PETSc.IntType], n: int, cols: npt.NDArray[PETSc.IntType], data: npt.NDArray[PETSc.ScalarTYpe], mode: int): MatSetValuesLocal(A, m, rows.ctypes, n, cols.ctypes, data.ctypes, mode) - MatSetValuesBlockedLocal = <_FuncPtr object>
- See PETSc MatSetValuesBlockedLocal documentation. 
 - MatSetValuesLocal = <_FuncPtr object>
- See PETSc MatSetValuesLocal documentation. 
 
- 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\). 
 
- class dolfinx.fem.petsc.numba_utils[source]
- Bases: - object- Utility attributes for working with Numba and PETSc. - These attributes are convenience functions for calling PETSc C functions from within Numba functions. - Note - Numba must be available to use these utilities. - Examples - A typical use of these utility functions is: - import numpy as np import numpy.typing as npt def set_vals(A: int, m: int, rows: npt.NDArray[PETSc.IntType], n: int, cols: npt.NDArray[PETSc.IntType], data: npt.NDArray[PETSc.ScalarTYpe], mode: int): MatSetValuesLocal(A, m, rows.ctypes, n, cols.ctypes, data.ctypes, mode) - MatSetValuesBlockedLocal = ExternalFunction(MatSetValuesBlockedLocal)
- See PETSc MatSetValuesBlockedLocal documentation. 
 - MatSetValuesLocal = ExternalFunction(MatSetValuesLocal)
- See PETSc MatSetValuesLocal documentation. 
 
- dolfinx.fem.petsc.set_bc(b: Vec, bcs: Sequence[DirichletBC] | Sequence[Sequence[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- gis 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 - bis nested or blocked,- bcsis a 2D array and- bcs[i]are the boundary conditions to apply to block/nest- i. Otherwise- bcsshould be a sequence of- DirichletBCs. For block/nest problems,- dolfinx.fem.bcs_by_block()can be used to prepare the 2D array of- DirichletBCobjects.
- x0 – Vector used in the value that constrained entries are set to. If - None,- x0is treated as zero.
- alpha – Scalar value used in the value that constrained entries are set to.