hifir4py.HIF

class HIF(A=None, **kw)[source]

The HIF preconditioner object

This class implements the flexible interface of the HIF preconditioner in Python, with supports of both mixed-precision computation and complex arithmetic.

Examples

One can simply create an empty preconditioner, and initialize it later.

>>> from hifir4py import *
>>> hif = HIF()

Alternatively, one can create an instance and factorize it.

>>> from scipy.sparse import rand
>>> A = rand(10, 10, 0.5)
>>> hif = HIF(A)

The latter usages shares the same interface as factorize does.

property A

User-input CRS (CSR) matrix

Type:csr_matrix
property M_

Underlying C++ HIF object

Warning

Access this only if you know what you are doing!

property S

User-input CRS (CSR) sparsifier

Type:csr_matrix
__init__(A=None, **kw)[source]

Create a HIF preconditioner

One can construct an empty preconditioner, i.e.,

>>> M = HIF()

Alternatively, one can construct and factorize at the same time, and see factorize for the interface.

__repr__() str[source]

Representation of a HIF instance

__str__() str[source]

Return str(self).

apply(b: numpy.ndarray, **kw) numpy.ndarray[source]

Apply the preconditioner with a given operation (op)

This function is to apply the preconditioner in the following four different modes:

  1. multilevel triangular solve (“S”),
  2. transpose/Hermitian multilevel triangular solve (“SH” or “ST”),
  3. multilevel matrix-vector multiplication (“M”), and
  4. tranpose/Hermitian multilevel matrix-vector multiplication (“MT” or “MH”).
Parameters:
  • b (ndarray) – Input RHS vector
  • op (str, optional) – Operation option, in {“S”, “SH”, “ST”, “M”, “MH”, “MT”} and case insensitive; note that “SH” and “ST” are the same
  • x (ndarray, optional) – Output result, can be passed in as workspace
  • betas (float or list, optional) – The residual bound for iterative refinement. If a scalar is passed in, then it is assumed to be lower bound (\(\beta_L\)) and the upper bound (\(\beta_U\)) is set to be maximum value.
  • nirs (int, optional) – Number of iterative refinement steps (default 1)
  • rank (int, optional) – Numerical rank used in final Schur complement
Returns:

  • x (ndarray) – Computed solution vector
  • iters (int, optional) – If iterative refinement is enabled in triangular solve and residual bounds betas is passed in, then this indicates the actual refinement iterations.
  • flag (int, optional) – If iterative refinement is enabled in triangular solve and residual bounds betas is passed in, then this indicates the status of iterative refinement. If flag==0, then the IR process converged; if flag>0, then it diverged; otherwise, it reached maximum iterative refinement limit (nirs).

Examples

Given an instance of HIF, say hif. The following computes the standard triangular solve in preconditioners

>>> x = hif.apply(b)

To enable transpose/Hermitian, one can do

>>> x = hif.apply(b, op="ST")

For triangular solve, the following enables two-step IR

>>> x = hif.apply(b, nirs=2)

Finally, the following demonstrates how to apply matrix-vector product operation

>>> x = hif.apply(b, op="M")
>>> x = hif.apply(b, op="MT")  # tranpose/Hermitian

See also

factorize, schur_rank

clear() None[source]

Clear the internal memory for the preconditioner

empty() bool[source]

Check emptyness

factorize(A, **kw) None[source]

Factorize a HIF preconditioner

This function is the core in HIF to (re)factorize a HIF preconditioner given input a matrix or sparsifier.

Parameters:
  • A (csr_matrix) – Input CRS matrix
  • S (csr_matrix or None, optional) – Optional sparsifier input (on which we will compute HIF)
  • is_mixed (bool, optional) – Whether or not using mixed-precision (using single)
  • params (Params, optional) – Control parameters, using default values if not provided

Examples

>>> from scipy.sparse import rand
>>> from hifir4py import *
>>> A = rand(10, 10, 0.5)
>>> hif = HIF()
>>> hif.factorize(A)

Notes

Besides passing in a Params object for control parameters, one can also directly using key-value pairs for parameters while calling the factorize function.

See also

apply

index_size() int[source]

Check the integer byte size used in underlying C++ HIF

Raises:AttributeError – This is raised if the underlying C++ HIF attribute is missing
is_complex() bool[source]

Check if the underlying C++ HIF is complex

Raises:AttributeError – This is raised if the underlying C++ HIF attribute is missing

See also

is_mixed

is_mixed() bool[source]

Check if the underlying HIF is mixed-precision

Raises:AttributeError – This is raised if the underlying C++ HIF attribute is missing

See also

is_complex

property levels: int

Number of levels

Type:int
property ncols: int

Number of columns

Type:int
property nnz: int

Total number of nonzeros in the preconditioner

Type:int
property nnz_ef: int

Total number of nonzeros in the E and F off diagonal blocks

Type:int
property nnz_ldu: int

Total number of nonzeros in the L, D, and U factors

Type:int
property nrows: int

Number of rows

Type:int
property rank: int

Numerical rank of the preconditioner

Type:int
refactorize(S, **kw) None[source]

Refactorize a sparsifier

Note

This function is similar to factorize, but differs in that this function doesn’t update the A matrix.

Parameters:
  • S (csr_matrix) – Sparsifier input (on which we will compute HIF)
  • is_mixed (bool, optional) – Whether or not using mixed-precision (using single)
  • params (Params, optional) – Control parameters, using default values if not provided

See also

update

Examples

>>> from scipy.sparse import rand
>>> from hifir4py import *
>>> hif = HIF(rand(10, 10, 0.5))
>>> A1 = hif.A
>>> A2 = rand(10, 10, 0.5)
>>> hif.refactorize(A2)
>>> assert A1 is hif.A

The following example illustrates how to update IR operator A and factorization asynchronously. In particular, we refactorize every 10 iterations, but update IR operator for every step.

>>> import numpy as np
>>> from scipy.sparse import rand
>>> from hifir4py import *
>>> hif = HIF(rand(10, 10, 0.5)) # initial HIF
>>> b = np.random.rand(10)
>>> for i in range(100):
>>>     x = hif.apply(b, nirs=2)  # 2-iter refinement with hif.A
>>>     A = rand(10, 10, 0.5)
>>>     hif.update(A)  # update A, or equiv as hif.A = A
>>>     if i % 10 == 0:
>>>         hif.refactorize(A)  # refactorization
property schur_rank: int

Numerical rank of the final Schur complement

Type:int
property schur_size: int

Size of the final Schur complement

Type:int
property shape: Tuple[int, int]

2-tuple of the preconditioner shape, i.e., (nrows, ncols)

Type:tuple
to_scipy()[source]

Compute a SciPy LinearOperator based on HIF

Returns:Return a linear operator in SciPy so that HIF can be used in its built-in KSP solvers.
Return type:LinearOperator
update(A) None[source]

Update the A matrix used in IR

Note

This function does not perform factorization on A.