Previous topic

mmf.math.linalg.cholesky.se99

This Page

mmf.math.linalg.pardiso

PardisoError Error with Pardiso solver.
Pardiso(*varargin, **kwargs) Encapsulate the Pardiso solver.

Inheritance diagram for mmf.math.linalg.pardiso:

Inheritance diagram of mmf.math.linalg.pardiso

Module with an interface to the Pardiso sparse solver.

>>> p = Pardiso()
>>> indptr = np.array([0, 4, 7, 9, 11, 14, 16, 17, 18])
>>> indices = np.array([0,    2,       5, 6,
...                        1, 2,    4,
...                           2,             7,
...                              3,       6,
...                                 4, 5, 6,
...                                    5,    7,
...                                       6,
...                                          7])
>>> data = np.array([7.0,      1.0,           2.0, 7.0,
...                      -4.0, 8.0,           2.0,
...                            1.0,                     5.0,
...                                 7.0,           9.0,
...                                      5.0, 1.0, 5.0,
...                                           0.0,      5.0,
...                                               11.0,
...                                                     5.0])
>>> A = sp.sparse.csr_matrix((data, indices, indptr))
>>> b = np.ones(A.shape[0], dtype=A.dtype)
>>> p.solve(A=A, b=b, symmetric=True)
...                              
array([-0.0427..., -0.0158... ,  0.1105..., -0.1204...,  0.0260...,
       -0.1224...,  0.2047...,  0.2118...])

An asymmetric example:

>>> indptr = np.array([0, 4, 7, 9, 11, 12, 15, 17, 20])
>>> indices = np.array([0,    2,       5, 6, 
...                        1, 2,    4,
...                           2,             7,
...                              3,       6,
...                        1,
...                           2,       5,    7,
...                        1,             6,
...                           2,          6, 7])
>>> data = np.array([7.0,      1.0,           2.0, 7.0, 
...                      -4.0, 8.0,      2.0,
...                            1.0,                     5.0,
...                                 7.0,           9.0,
...                      -4.0,
...                            7.0,           3.0,      8.0,
...                       1.0,                    11.0,
...                           -3.0,                2.0, 5.0])
>>> p.solve(A=A, b=b, symmetric=False)
...                              
array([ 0. , -0.213...,  0. ,  0.0259...,  0.07272...,
        0.1818...,  0.0909...,  0.2 ])
class mmf.math.linalg.pardiso.PardisoError[source]

Bases: exceptions.Exception

Error with Pardiso solver.

__init__()

x.__init__(...) initializes x; see help(type(x)) for signature

class mmf.math.linalg.pardiso.Pardiso(*varargin, **kwargs)[source]

Bases: mmf.objects.StateVars

Encapsulate the Pardiso solver.

Pardiso(solver=0,
        threads=None,
        verbose=False)

Notes

The solver expects the input for matrices in the scipy.sparse.csr_matrix format. If the matrices are symmetric, only the upper triangular portion should be provided.

Attributes

State Variables:  
solver

This scalar value defines the solver method that the user would like to use:

0 : Sparse direct solver 1 : Multi-recursive iterative solver

threads
Number of threads to use. If None, then look at the environment variable OMP_NUM_THREADS or default to 1.
verbose
If True, then pardiso will print out a buch of statistics.

Methods

all_items() Return a list [(name, obj)] of (name, object) pairs containing all items available.
archive(name) Return a string representation that can be executed to restore the object.
archive_1([env]) Return (rep, args, imports).
initialize(**kwargs) Calls __init__() passing all assigned attributes.
items([archive]) Return a list [(name, obj)] of (name, object) pairs where the
mtype(A[, structurally_symmetric, ...]) Return (mtype, A) for the matrix A.
resume() Resume calls to __init__() from __setattr__(),
solve(A, b, **kwargs) Return solution x to A*x = b.
suspend() Suspend calls to __init__() from
__init__(*varargin, **kwargs)[source]
mtype(A, structurally_symmetric=False, symmetric=False, hermitian=False, positive_definite=False, auto=False)[source]

Return (mtype, A) for the matrix A.

Parameters :

A : sparse array

The matrix. This is used for A.dtype and if auto is True.

structurally_symmetric : bool

Set to True if the sparsity pattern is symmetric.

symmetric : bool

Set to True if the matrix is symmetric. If so, only the upper triangular portion will be returned.

hermitian : bool

Set to True if the matrix is complex and Hermitian. If so, only the upper triangular portion will be returned.

positive_definite : bool

Set to True if the matrix is real and symmetric positive definite or complex and Hermitian positive definite. If so, only the upper triangular portion will be returned.

auto : bool

If True, Attempt to figure out the best structure for the matrix. Will assume not positive definite.

Returns :

mtype : int

Matrix type for pardiso library.

A : csr_matrix

Array packed in appropriate scipy.sparse.csr_matrix format. Only upper triangular portion is non-zero if Hermitian or symmetric.

solve(A, b, **kwargs)[source]

Return solution x to A*x = b.

Parameters :

A : array

Matrix. Will be converted to

scipy.sparse.csr_matrix so provide in this format if possible.

b : array

Right-hand side(s). If this has two-dimensions, then b.shape[1] equations will be solved simultaneously.

**kwargs : :

See also the parameters of mtype() for additional keyword arguments allowed for specifying the type of matrix.

Returns :

——- :

x : array

Solution to A*x = b.