dolfinx.la.superlu_dist

SuperLU_DIST linear solver support.

This module provides support for parallel solution of linear systems assembled into dolfinx.la.MatrixCSR using SuperLU_DIST.

Note

Users with advanced linear solver requirements should use PETSc/petsc4py.

Functions

superlu_dist_matrix(A)

Create a SuperLU_DIST matrix.

superlu_dist_solver(A)

Create a SuperLU_DIST linear solver.

Classes

SuperLUDistMatrix(matrix)

SuperLU_DIST matrix.

SuperLUDistSolver(solver)

SuperLU_DIST solver.

class dolfinx.la.superlu_dist.SuperLUDistMatrix(matrix)[source]

Bases: object

SuperLU_DIST matrix.

Create a SuperLU_DIST matrix.

Parameters:

matrix – C++ SuperLUDistMatrix object.

Note

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

property dtype: DTypeLike

Dtype of matrix values.

class dolfinx.la.superlu_dist.SuperLUDistSolver(solver)[source]

Bases: object

SuperLU_DIST solver.

Create a SuperLU_DIST solver.

Parameters:

solver – C++ SuperLUDistSolver object.

Note

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

set_A(A: SuperLUDistMatrix)[source]

Set assembled left-hand side matrix.

For advanced use with SuperLU_DIST option Factor allowing use of previously computed permutations when solving with new matrix A.

Parameters:

A – Assembled left-hand side matrix \(A\).

set_option(name: str, value: str)[source]

Set SuperLU_DIST option for solve.

See SuperLU_DIST User’s Guide for option names and values.

Examples

solver.set_option(“SymmetricMode”, “YES”) solver.set_option(“Trans”, “NOTRANS”)

Parameters:
  • name – Option name.

  • value – Option value.

solve(b: Vector, u: Vector) int[source]

Solve linear system \(Au = b\).

Note

Vectors must have size and parallel layout compatible with A.

Note

The caller must check the return integer for success (== 0).

Note

The caller must u.scatter_forward() after the solve.

Note

The values of A are modified in-place during the solve.

Note

To solve with successive right-hand sides b the caller must solver.set_option("Factor", "FACTORED") after the first solve.

Parameters:
  • b – Right-hand side vector \(b\).

  • u – Solution vector \(u\), overwritten during solve.

Returns:

SuperLU_DIST return integer from p*gssvx routine.

dolfinx.la.superlu_dist.superlu_dist_matrix(A: MatrixCSR) SuperLUDistMatrix[source]

Create a SuperLU_DIST matrix.

Deep copies all required data from A.

Parameters:

A – Assembled matrix.

Returns:

A SuperLU_DIST matrix.

dolfinx.la.superlu_dist.superlu_dist_solver(A: SuperLUDistMatrix) SuperLUDistSolver[source]

Create a SuperLU_DIST linear solver.

Solve linear system \(Au = b\) via LU decomposition.

The SuperLU_DIST solver has options set to upstream defaults, except PrintStat (verbose solver output) set to NO.

Parameters:

A – Assembled left-hand side matrix \(A\).

Returns:

A SuperLU_DIST solver.