1
2PRIMME: PReconditioned Iterative MultiMethod Eigensolver
3========================================================
4
5`primme` is a Python interface to PRIMME_, a high-performance library for computing a few eigenvalues/eigenvectors, and singular values/vectors.
6PRIMME is especially optimized for large, difficult problems.
7Real symmetric and complex Hermitian problems, standard :math:`A x = \lambda x` and generalized :math:`A x = \lambda B x`, are supported.
8It can find largest, smallest, or interior singular/eigenvalues, and can use preconditioning to accelerate convergence.
9
10The main contributors to PRIMME are James R. McCombs, Eloy Romero Alcalde, Andreas Stathopoulos and Lingfei Wu.
11
12Install
13-------
14
15You can install the latest version with `pip`::
16
17    pip install numpy   # if numpy is not installed yet
18    pip install scipy   # if scipy is not installed yet
19    pip install future  # if using python 2
20    conda install mkl-devel # if using Anaconda Python distribution
21    pip install primme
22
23Optionally for building the development version do::
24
25    git clone https://github.com/primme/primme
26    cd primme
27    make python_install
28
29Usage
30-----
31
32The following examples compute a few eigenvalues and eigenvectors from a real symmetric matrix::
33
34    >>> import Primme, scipy.sparse
35    >>> A = scipy.sparse.spdiags(range(100), [0], 100, 100) # sparse diag. matrix
36    >>> evals, evecs = Primme.eigsh(A, 3, tol=1e-6, which='LA')
37    >>> evals # the three largest eigenvalues of A
38    array([ 99.,  98.,  97.])
39
40    >>> new_evals, new_evecs = Primme.eigsh(A, 3, tol=1e-6, which='LA', ortho=evecs)
41    >>> new_evals # the next three largest eigenvalues
42    array([ 96.,  95.,  94.])
43
44    >>> evals, evecs = primme.eigsh(A, 3, tol=1e-6, which=50.1)
45    >>> evals # the three closest eigenvalues to 50.1
46    array([ 50.,  51.,  49.])
47
48
49The following examples compute a few eigenvalues and eigenvectors from a generalized Hermitian problem, without factorizing or inverting :math:`B`::
50
51    >>> import Primme, scipy.sparse
52    >>> A = scipy.sparse.spdiags(range(100), [0], 100, 100) # sparse diag. matrix
53    >>> M = scipy.sparse.spdiags(np.asarray(range(99,-1,-1)), [0], 100, 100)
54    >>> evals, evecs = primme.eigsh(A, 3, M=M, tol=1e-6, which='SA')
55    >>> evals
56    array([1.0035e-07, 1.0204e-02, 2.0618e-02])
57
58The following examples compute a few singular values and vectors::
59
60    >>> import Primme, scipy.sparse
61    >>> A = scipy.sparse.spdiags(range(1, 11), [0], 100, 10) # sparse diag. rect. matrix
62    >>> svecs_left, svals, svecs_right = Primme.svds(A, 3, tol=1e-6, which='SM')
63    >>> svals # the three smallest singular values of A
64    array([ 1.,  2.,  3.])
65
66    >>> A = scipy.sparse.rand(10000, 100, random_state=10)
67    >>> prec = scipy.sparse.spdiags(np.reciprocal(A.multiply(A).sum(axis=0)),
68    ...           [0], 100, 100) # square diag. preconditioner
69    >>> svecs_left, svals, svecs_right = Primme.svds(A, 3, which=6.0, tol=1e-6,
70    ...           precAHA=prec)
71    >>> ["%.5f" % x for x in svals.flat] # the three closest singular values of A to 0.5
72    ['5.99871', '5.99057', '6.01065']
73
74Further examples_.
75
76Documentation of eigsh_ and svds_.
77
78Citing this code
79----------------
80
81Please cite (bibtex_):
82
83* A. Stathopoulos and J. R. McCombs *PRIMME: PReconditioned Iterative
84  MultiMethod Eigensolver: Methods and software description*, ACM
85  Transaction on Mathematical Software Vol. 37, No. 2, (2010),
86  21:1-21:30.
87
88* L. Wu, E. Romero and A. Stathopoulos, *PRIMME_SVDS: A High-Performance
89  Preconditioned SVD Solver for Accurate Large-Scale Computations*,
90  J. Sci. Comput., Vol. 39, No. 5, (2017), S248--S271.
91
92License Information
93-------------------
94
95PRIMME and this interface is licensed under the 3-clause license BSD.
96
97Contact Information
98-------------------
99
100For reporting bugs or questions about functionality contact `Andreas Stathopoulos`_ by
101email, `andreas` at `cs.wm.edu`. See further information in
102the webpage http://www.cs.wm.edu/~andreas/software.
103
104.. _PRIMME: https://github.com/primme/primme
105.. _`Andreas Stathopoulos`: http://www.cs.wm.edu/~andreas/software
106.. _`github`: https://github.com/primme/primme
107.. _`doc`: http://www.cs.wm.edu/~andreas/software/doc/readme.html
108.. _PETSc : http://www.mcs.anl.gov/petsc/
109.. _`bibtex`: https://raw.githubusercontent.com/primme/primme/master/doc/primme.bib
110.. _eigsh: http://www.cs.wm.edu/~andreas/software/doc/pyeigsh.html
111.. _svds: http://www.cs.wm.edu/~andreas/software/doc/pysvds.html
112.. _examples: https://github.com/primme/primme/blob/master/Python/examples.py
113