dolfinx.la

Linear algebra functionality

Functions

orthonormalize(basis)

Orthogonalise set of PETSc vectors in-place.

is_orthonormal(basis[, eps])

Check that list of PETSc vectors are orthonormal.

matrix_csr(sp[, block_mode, dtype])

Create a distributed sparse matrix.

vector(map[, bs, dtype])

Create a distributed vector.

create_petsc_vector(map, bs)

Create a distributed PETSc vector.

Classes

MatrixCSR(A)

A distributed sparse matrix that uses compressed sparse row storage.

Norm

InsertMode

Vector(x)

A distributed vector object.

class dolfinx.la.InsertMode

Bases: object

add = dolfinx.cpp.la.InsertMode.add
insert = dolfinx.cpp.la.InsertMode.insert
class dolfinx.la.MatrixCSR(A: MatrixCSR_float32 | MatrixCSR_float64 | MatrixCSR_complex64 | MatrixCSR_complex128)[source]

Bases: object

A distributed sparse matrix that uses compressed sparse row storage.

Note

Objects of this type should be created using matrix_csr() and not created using this initialiser.

Parameters:

A – The C++/nanobind matrix object.

add(x: ndarray[Any, dtype[floating]], rows: ndarray[Any, dtype[int32]], cols: ndarray[Any, dtype[int32]], bs: int = 1) None[source]

Add a block of values in the matrix.

property block_size

Block sizes for the matrix.

property data: ndarray[Any, dtype[floating]]

Underlying matrix entry data.

index_map(i: int) IndexMap[source]

Index map for row/column.

Parameters:

i – 0 for row map, 1 for column map.

property indices: ndarray[Any, dtype[int32]]

Local column indices.

property indptr: ndarray[Any, dtype[int64]]

Local row pointers.

scatter_reverse() None[source]

Scatter and accumulate ghost values.

set(x: ndarray[Any, dtype[floating]], rows: ndarray[Any, dtype[int32]], cols: ndarray[Any, dtype[int32]], bs: int = 1) None[source]

Set a block of values in the matrix.

set_value(x: floating) None[source]

Set all non-zero entries to a value.

Parameters:

x – The value to set all non-zero entries to.

squared_norm() floating[source]

Compute the squared Frobenius norm.

Note

This operation is collective and requires communication.

to_dense() ndarray[Any, dtype[floating]][source]

Copy to a dense 2D array.

Note

Typically used for debugging.

to_scipy(ghosted=False)[source]

Convert to a SciPy CSR/BSR matrix. Data is shared.

Note

SciPy must be available.

Parameters:

ghosted – If True rows that are ghosted in parallel are included in the returned SciPy matrix, otherwise ghost rows are not included.

Returns:

SciPy compressed sparse row (both block sizes equal to one) or a SciPy block compressed sparse row matrix.

class dolfinx.la.Norm

Bases: object

frobenius = dolfinx.cpp.la.Norm.frobenius
l1 = dolfinx.cpp.la.Norm.l1
l2 = dolfinx.cpp.la.Norm.l2
linf = dolfinx.cpp.la.Norm.linf
class dolfinx.la.Vector(x: Vector_float32 | Vector_float64 | Vector_complex64 | Vector_complex128 | Vector_int8 | Vector_int32 | Vector_int64)[source]

Bases: object

A distributed vector object.

Parameters:

x – C++ Vector object.

Note

This initialiser is intended for internal library use only. User code should call vector() to create a vector object.

property array: ndarray

Local representation of the vector.

property block_size: int

Block size for the vector.

property index_map: IndexMap

Index map that describes size and parallel distribution.

property petsc_vec

PETSc vector holding the entries of the vector.

Upon first call, this function creates a PETSc Vec object that wraps the degree-of-freedom data. The Vec object is cached and the cached Vec is returned upon subsequent calls.

scatter_forward() None[source]

Update ghost entries.

scatter_reverse(mode: InsertMode) None[source]

Scatter ghost entries to owner.

Parameters:

mode – Control how scattered values are set/accumulated by owner.

dolfinx.la.create_petsc_vector(map, bs: int)[source]

Create a distributed PETSc vector.

Parameters:
  • map – Index map that describes the size and parallel layout of the vector to create.

  • bs – Block size of the vector.

Returns:

PETSc Vec object.

dolfinx.la.is_orthonormal(basis, eps: float = 1e-12) bool[source]

Check that list of PETSc vectors are orthonormal.

dolfinx.la.matrix_csr(sp: ~dolfinx.cpp.la.SparsityPattern, block_mode=dolfinx.cpp.la.BlockMode.compact, dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>) MatrixCSR[source]

Create a distributed sparse matrix.

The matrix uses compressed sparse row storage.

Parameters:
  • sp – The sparsity pattern that defines the nonzero structure of the matrix the parallel distribution of the matrix.

  • dtype – The scalar type.

Returns:

A sparse matrix.

dolfinx.la.orthonormalize(basis)[source]

Orthogonalise set of PETSc vectors in-place.

dolfinx.la.vector(map, bs=1, dtype: ~numpy.dtype[~typing.Any] | None | type[~typing.Any] | ~numpy._typing._dtype_like._SupportsDType[~numpy.dtype[~typing.Any]] | str | tuple[~typing.Any, int] | tuple[~typing.Any, ~typing.SupportsIndex | ~collections.abc.Sequence[~typing.SupportsIndex]] | list[~typing.Any] | ~numpy._typing._dtype_like._DTypeDict | tuple[~typing.Any, ~typing.Any] = <class 'numpy.float64'>) Vector[source]

Create a distributed vector.

Parameters:
  • map – Index map the describes the size and distribution of the vector.

  • bs – Block size.

  • dtype – The scalar type.

Returns:

A distributed vector.