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()

Assign x0 values to a PETSc vector x1.

create_vector(maps[, kind])

Create a PETSc vector from a sequence of maps and blocksizes.

create_vector_wrap(x)

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 x0 values to a PETSc vector x1.

Values in x0, which is possibly a stacked collection of arrays, are assigned x1. When x0 holds a sequence of n` arrays and x1 has type NEST, the assignment is:

      [x0[0]]
x1 =  [x0[1]]
      [.....]
      [x0[n-1]]

When x0 holds a sequence of n arrays and x1 does not have type NEST, 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:

  1. If maps=[(im_0, bs_0), ..., (im_n, bs_n)] is a sequence of indexmaps and blocksizes and kind is None``or is ``PETSc.Vec.Type.MPI, a ghosted PETSc vector whith block structure described by (im_i, bs_i) 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 (im_i, bs_i) and b_ig are the ‘unowned’ (ghost) entries.

    If more than one tuple is supplied, the returned vector has an attribute _blocks that holds the local offsets into b for 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 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]]
    
  2. If V=[(im_0, bs_0), ..., (im_n, bs_n)] is a sequence of function space and kind is PETSc.Vec.Type.NEST, a PETSc nested vector (a ‘nest’ of ghosted PETSc vectors) is created.

Parameters:
  • map – Sequence of tuples of IndexMap and 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.

dolfinx.la.petsc.create_vector_wrap(x: Vector) Vec[source]

Wrap a distributed DOLFINx vector as a PETSc vector.

Parameters:

x – The vector to wrap as a PETSc vector.

Returns:

A PETSc vector that shares data with x.