cudolfinx.assemble
GPU-accelerated assembly routines.
Functions
Create PETSc Vector on device. |
Classes
Class for assembly on the GPU. |
- class cudolfinx.assemble.CUDAAssembler[source]
Bases:
objectClass for assembly on the GPU.
Initialize the assembler.
- apply_lifting(b: CUDAVector, a: list[CUDAForm], bcs: list[list[DirichletBC] | CUDADirichletBC], x0: list[CUDAVector] | None = None, scale: float = 1.0, coeffs: list[list[Function]] | None = None)[source]
GPU equivalent of apply_lifting.
- Parameters:
b – CUDAVector to modify
a – list of forms to lift
bcs – list of boundary condition lists
x0 – optional list of shift vectors for lifting
scale – scale of lifting
coeffs – coefficients to (re-)pack
- apply_lifting_block(b: CUDAVector, a: BlockCUDAForm, bcs: list[DirichletBC] | Any, x0: CUDAVector | None = None, alpha: float = 1.0, set_bcs=True)[source]
Apply lifting with a single block form.
- assemble_matrix(a: CUDAForm, mat: CUDAMatrix | None = None, bcs: list[DirichletBC] | CUDADirichletBC | None = [], diagonal: float = 1.0, constants: list | None = None, coeffs: list | None = None, zero: bool = True)[source]
Assemble bilinear form into a matrix on the GPU.
- Parameters:
a – The bilinear form to assemble.
mat – Matrix in which to store assembled values.
bcs – Boundary conditions that affect the assembled matrix. Degrees-of-freedom constrained by a boundary condition will have their rows/columns zeroed and the value
diagonalset on on the matrix diagonal.diagonal – Value to set on the diagonal.
constants – Constants that appear in the form. If not provided, any required coefficients will be computed.
coeffs – Optional list of form coefficients to repack. If not specified, all will be repacked. If an empty list is passed, no repacking will be performed.
zero – Whether to zero the matrix prior to assembly
- Returns:
Matrix representation of the bilinear form
a.
Note
The returned matrix is not finalised, i.e. ghost values are not accumulated.
- assemble_matrix_block(a: BlockCUDAForm, mat: CUDAMatrix | None = None, bcs: list[Any] | None = None, coeffs: list[list[Any]] | None = None, constants: list[list[Any]] | None = None, zero: bool = True)[source]
Assemble block form into a matrix on the GPU.
- assemble_scalar(b: CUDAForm, constants=None, coeffs=None)[source]
Assemble scalar integral on GPU.
- Parameters:
b – the functional to use for assembly
constants – Form constants
coeffs – Optional list of form coefficients to repack. If not specified, all will be repacked. If an empty list is passed, no repacking will be performed.
- assemble_vector(b: CUDAForm, vec: CUDAVector | None = None, constants=None, coeffs=None, zero=True)[source]
Assemble linear form into vector on GPU.
- Parameters:
b – the linear form to use for assembly
vec – the vector to assemble into. Created if it doesn’t exist
constants – Form constants
coeffs – Optional list of form coefficients to repack. If not specified, all will be repacked. If an empty list is passed, no repacking will be performed.
zero – Whether to zero vector entries before assembly (defaults to true)
- assemble_vector_block(b: BlockCUDAForm, vec: CUDAVector | None = None, constants=None, coeffs=None, zero=True)[source]
Assemble block linear form into vector on GPU.
- Parameters:
b – the block linear form to use for assembly
vec – the vector to assemble into. Created if not provided
constants – List of constants for each sub form
coeffs –
- Optional list of form coefficients to repack.
If not specified, all will be repacked.
If an empty list is passed, no repacking will be performed.
zero – Whether to zero vector entries before assembly (defaults to true)
- create_matrix(a: CUDAForm) CUDAMatrix[source]
Create a CUDAMatrix from a given form.
- create_matrix_block(a: BlockCUDAForm) CUDAMatrix[source]
Create a block matrix from a block form.
- create_vector(b: CUDAForm) CUDAVector[source]
Create a CUDAVector from a given form.
- create_vector_block(b: BlockCUDAForm) CUDAVector[source]
Create a CUDAVector from a given form.
- pack_bcs(bcs: list[DirichletBC]) CUDADirichletBC[source]
Pack boundary conditions into a single object for use in assembly.
The returned object is of type CUDADirichletBC and can be used in place of a list of regular DirichletBCs. This is more efficient when performing multiple operations with the same list of boundary conditions, or when boundary condition values need to change over time.
- pack_coefficients(a: CUDAForm, coefficients: list[Function] | None = None)[source]
Pack coefficients on device.
- set_bc(b: CUDAVector, bcs: list[DirichletBC] | CUDADirichletBC, V: FunctionSpace, x0: CUDAVector | None = None, scale: float = 1.0)[source]
Set boundary conditions on device.
- Parameters:
b – vector to modify
bcs – collection of bcs
V – FunctionSpace to which bcs apply
x0 – optional shift vector
scale – scaling factor