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 a PETSc vector that is compaible with a linear form. |
|
Create a PETSc vector (blocked) that is compaible with a list of linear forms. |
|
Create a PETSc nested vector ( |
|
|
Create a PETSc matrix that is compatible with a bilinear form. |
Create a PETSc matrix that is compatible with a rectangular array of bilinear forms. |
|
Create a PETSc matrix ( |
|
Assemble linear form into a new PETSc vector. |
|
Assemble linear forms into a new nested PETSc ( |
|
Assemble linear forms into a monolithic vector. |
|
Assemble bilinear form into a matrix. |
|
Create a nested matrix and assembled bilinear forms into the matrix. |
|
Assemble bilinear forms into a blocked matrix. |
|
|
Apply the function |
|
Apply the function |
|
Apply the function |
|
Apply the function |
Classes
|
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. |
|
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: Form, L: Form, bcs: List[DirichletBC] = [], u: Optional[Function] = None, petsc_options: Optional[dict] = None, form_compiler_options: Optional[dict] = None, jit_options: Optional[dict] = None)[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 – 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 b: Vec
Right-hand side vector
- property solver: KSP
Linear solver object
- class dolfinx.fem.petsc.NonlinearProblem(F: Form, u: Function, bcs: List[DirichletBC] = [], J: Optional[Form] = None, form_compiler_options: Optional[dict] = None, jit_options: Optional[dict] = None)[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_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(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
- dolfinx.fem.petsc.apply_lifting(b: Vec, a: List[Form], bcs: List[List[DirichletBC]], 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[Form]], bcs: List[DirichletBC], 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[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 matrix diagonal for Dirichlet boundary condition constrained degrees-of-freedom.
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 assembled 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 matrix diagonal for Dirichlet boundary condition constrained degrees-of-freedom.
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: 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(b: Vec, L: List[Form], a: List[List[Form]], bcs: List[DirichletBC] = [], x0: Optional[Vec] = None, scale: float = 1.0, 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.
- 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.
- 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 a rectangular array of bilinear forms.- Parameters:
a – Rectangular array of bilinear forms.
- Returns:
A PETSc matrix (
MatNest
) that is compatible witha
.
- dolfinx.fem.petsc.create_vector(L: Form) 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[Form]) 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[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 withL
.
- dolfinx.fem.petsc.set_bc(b: Vec, bcs: List[DirichletBC], 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[DirichletBC]], 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.