dolfinx.la.petsc
Functions for working with PETSc linear algebra objects.
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
|
Assign |
|
Create a PETSc vector from a sequence of maps and blocksizes. |
Wrap a distributed DOLFINx vector as a PETSc vector. |
- dolfinx.la.petsc.assign(x0: ndarray[tuple[int, ...], dtype[inexact]] | Sequence[ndarray[tuple[int, ...], dtype[inexact]]], x1: Vec)[source]
- dolfinx.la.petsc.assign(x0: Vec, x1: ndarray[tuple[int, ...], dtype[inexact]] | Sequence[ndarray[tuple[int, ...], dtype[inexact]]])
Assign
x0values to a PETSc vectorx1.Values in
x0, which is possibly a stacked collection of arrays, are assignedx1. Whenx0holds a sequence ofn`arrays andx1has typeNEST, the assignment is:[x0[0]] x1 = [x0[1]] [.....] [x0[n-1]]
When
x0holds a sequence ofnarrays andx1does not have typeNEST, the assignment is:[x0_owned[0]] x1 = [.....] [x0_owned[n-1]] [x0_ghost[0]] [.....] [x0_ghost[n-1]]
- Parameters:
x0 – An array or list of arrays that will be assigned to
x1.x1 – Vector to assign values to.
- dolfinx.la.petsc.create_vector(maps: Sequence[tuple[IndexMap, int]], kind: str | None = None) Vec[source]
Create a PETSc vector from a sequence of maps and blocksizes.
Three cases are supported:
If
maps=[(im_0, bs_0), ..., (im_n, bs_n)]is a sequence of indexmaps and blocksizes andkindisNone``or is ``PETSc.Vec.Type.MPI, a ghosted PETSc vector whith block structure described by(im_i, bs_i)is created. The created vectorbis initialized such that on each MPI processb = [b_0, b_1, ..., b_n, b_0g, b_1g, ..., b_ng], whereb_iare the entries associated with the ‘owned’ degrees-of- freedom for(im_i, bs_i)andb_igare the ‘unowned’ (ghost) entries.If more than one tuple is supplied, the returned vector has an attribute
_blocksthat holds the local offsets intobfor the (i) owned and (ii) ghost entries for eachV[i]. It can be accessed byb.getAttr("_blocks"). The offsets can be used to get views intobfor 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
V=[(im_0, bs_0), ..., (im_n, bs_n)]is a sequence of function space andkindisPETSc.Vec.Type.NEST, a PETSc nested vector (a ‘nest’ of ghosted PETSc vectors) is created.
- Parameters:
map – Sequence of tuples of
IndexMapand the associated block size.kind – PETSc vector type (
VecType) to create.
- Returns:
A PETSc vector with the prescribed layout. The vector is not initialised to zero.