File: arpack.py

package info (click to toggle)
python-scipy 0.5.2-0.1
  • links: PTS
  • area: main
  • in suites: etch, etch-m68k
  • size: 33,888 kB
  • ctags: 44,231
  • sloc: ansic: 156,256; cpp: 90,347; python: 89,604; fortran: 73,083; sh: 1,318; objc: 424; makefile: 342
file content (383 lines) | stat: -rw-r--r-- 11,625 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
"""
arpack - Scipy module to find a few eigenvectors and eigenvalues of a matrix

Uses ARPACK: http://www.caam.rice.edu/software/ARPACK/

"""
__all___=['eigen','eigen_symmetric']

import _arpack 
import numpy as sb
import warnings

# inspired by iterative.py
# so inspired, in fact, that some of it was copied directly
try:
    False, True
except NameError:
    False, True = 0, 1

_type_conv = {'f':'s', 'd':'d', 'F':'c', 'D':'z'}

class get_matvec:
    methname = 'matvec'
    def __init__(self, obj, *args):
        self.obj = obj
        self.args = args
        if isinstance(obj, sb.matrix):
            self.callfunc = self.type1m
            return
        if isinstance(obj, sb.ndarray):
            self.callfunc = self.type1
            return
        meth = getattr(obj,self.methname,None)
        if not callable(meth):
            raise ValueError, "Object must be an array "\
                  "or have a callable %s attribute." % (self.methname,)

        self.obj = meth
        self.callfunc = self.type2

    def __call__(self, x):
        return self.callfunc(x)

    def type1(self, x):
        return sb.dot(self.obj, x)

    def type1m(self, x):
        return sb.dot(self.obj.A, x)

    def type2(self, x):
        return self.obj(x,*self.args)


def eigen(A,k=6,M=None,ncv=None,which='LM',
          maxiter=None,tol=0, return_eigenvectors=True):
    """ Return k eigenvalues and eigenvectors of the matrix A.

    Solves A * x[i] = w[i] * x[i], the standard eigenvalue problem for
    w[i] eigenvalues with corresponding eigenvectors x[i].

    Inputs:

    A --  A matrix, array or an object with matvec(x) method to perform
          the matrix vector product A * x.  The sparse matrix formats
          in scipy.sparse are appropriate for A.

    k -- The number of eigenvalue/eigenvectors desired

    M -- (Not implemented)
          A symmetric positive-definite matrix for the generalized
          eigenvalue problem A * x = w * M * x

    Outputs:

    w -- An array of k eigenvalues

    v -- An array of k eigenvectors, k[i] is the eigenvector corresponding
         to the eigenvector w[i]

    Optional Inputs:

    ncv -- Number of Lanczos vectors generated, ncv must be greater than k
           and is recommended to be ncv > 2*k

    which -- String specifying which eigenvectors to compute.
             Compute the k eigenvalues of:
             'LM' - largest magnitude.
             'SM' - smallest magnitude.
             'LR' - largest real part.
             'SR' - smallest real part.
             'LI' - largest imaginary part.
             'SI' - smallest imaginary part.

    maxiter --  Maximum number of Arnoldi update iterations allowed

    tol -- Relative accuracy for eigenvalues (stopping criterion)

    return_eigenvectors -- True|False, return eigenvectors 

    """
    A=sb.asanyarray(A)
    n,ny=A.shape
    try:
        n==ny
    except:
        raise AttributeError("matrix is not square")
    if M is not None:
        raise NotImplementedError("generalized eigenproblem not supported yet")
        
    # some defaults
    if ncv is None:
        ncv=2*k+1
    ncv=min(ncv,n)
    if maxiter==None:
        maxiter=n*10

    # guess type        
    resid = sb.zeros(n,'f')
    try:
        typ = A.dtype.char
    except AttributeError:
        typ = A.matvec(resid).dtype.char
    if typ not in 'fdFD':
        raise ValueError("matrix type must be 'f', 'd', 'F', or 'D'")

    # some sanity checks
    if k <= 0:
        raise ValueError("k must be positive, k=%d"%k)
    if k == n:
        raise ValueError("k must be less than rank(A), k=%d"%k)
    if maxiter <= 0:
        raise ValueError("maxiter must be positive, maxiter=%d"%maxiter)
    whiches=['LM','SM','LR','SR','LI','SI']
    if which not in whiches:
        raise ValueError("which must be one of %s"%' '.join(whiches))
    if ncv > n or ncv < k:
        raise ValueError("ncv must be k<=ncv<=n, ncv=%s"%ncv)

    # assign solver and postprocessor
    ltr = _type_conv[typ]
    eigsolver = _arpack.__dict__[ltr+'naupd']
    eigextract = _arpack.__dict__[ltr+'neupd']
    matvec = get_matvec(A)

    v = sb.zeros((n,ncv),typ) # holds Ritz vectors
    resid = sb.zeros(n,typ) # residual
    workd = sb.zeros(3*n,typ) # workspace
    workl = sb.zeros(3*ncv*ncv+6*ncv,typ) # workspace
    iparam = sb.zeros(11,'int') # problem parameters
    ipntr = sb.zeros(14,'int') # pointers into workspaces
    info = 0
    ido = 0

    if typ in 'FD':
        rwork = sb.zeros(ncv,typ.lower())

    # only supported mode is 1: Ax=lx
    ishfts = 1
    mode1 = 1 
    bmat = 'I'
    iparam[0] = ishfts
    iparam[2] = maxiter
    iparam[6] = mode1

    while True:
        if typ in 'fd':
            ido,resid,v,iparam,ipntr,info =\
                eigsolver(ido,bmat,which,k,tol,resid,v,iparam,ipntr,
                          workd,workl,info)
        else:
            ido,resid,v,iparam,ipntr,info =\
                eigsolver(ido,bmat,which,k,tol,resid,v,iparam,ipntr,
                          workd,workl,rwork,info)

        if (ido == -1 or ido == 1):
            # compute y = A * x
            xslice = slice(ipntr[0]-1, ipntr[0]-1+n)
            yslice = slice(ipntr[1]-1, ipntr[1]-1+n)
            workd[yslice]=matvec(workd[xslice])
        else: # done
            break

    if  info < -1 :
        raise RuntimeError("Error info=%d in arpack"%info)
        return None
    if info == -1:
        warnings.warn("Maximum number of iterations taken: %s"%iparam[2])
#    if iparam[3] != k:
#        warnings.warn("Only %s eigenvalues converged"%iparam[3])


    # now extract eigenvalues and (optionally) eigenvectors        
    rvec = return_eigenvectors
    ierr = 0
    howmny = 'A' # return all eigenvectors
    sselect = sb.zeros(ncv,'int') # unused 
    sigmai = 0.0 # no shifts, not implemented
    sigmar = 0.0 # no shifts, not implemented
    workev = sb.zeros(3*ncv,typ) 

    if typ in 'fd':
        dr=sb.zeros(k+1,typ)
        di=sb.zeros(k+1,typ)
        zr=sb.zeros((n,k+1),typ)
        dr,di,z,info=\
            eigextract(rvec,howmny,sselect,sigmar,sigmai,workev,
                   bmat,which,k,tol,resid,v,iparam,ipntr,
                   workd,workl,info)
        
        # make eigenvalues complex
        d=dr+1.0j*di
        # futz with the eigenvectors:
        # complex are stored as real,imaginary in consecutive columns
        z=zr.astype(typ.upper())
        for i in range(k): # fix c.c. pairs
            if di[i] > 0 :
                z[:,i]=zr[:,i]+1.0j*zr[:,i+1]
                z[:,i+1]=z[:,i].conjugate()
                     
    else:
        d,z,info =\
              eigextract(rvec,howmny,sselect,sigmar,workev,
                         bmat,which,k,tol,resid,v,iparam,ipntr,
                         workd,workl,rwork,ierr)

        

    if ierr != 0:
        raise RuntimeError("Error info=%d in arpack"%info)
        return None
    if return_eigenvectors:
        return d,z
    return d


def eigen_symmetric(A,k=6,M=None,ncv=None,which='LM',
                    maxiter=None,tol=0, return_eigenvectors=True):
    """ Return k eigenvalues and eigenvectors of the real symmetric matrix A.

    Solves A * x[i] = w[i] * x[i], the standard eigenvalue problem for
    w[i] eigenvalues with corresponding eigenvectors x[i].
    A must be real and symmetric.
    See eigen() for nonsymmetric or complex symmetric (Hermetian) matrices.

    Inputs:

    A --  A symmetric matrix, array or an object with matvec(x) method
          to perform the matrix vector product A * x.
          The sparse matrix formats in scipy.sparse are appropriate for A.

    k -- The number of eigenvalue/eigenvectors desired

    M -- (Not implemented)
          A symmetric positive-definite matrix for the generalized
          eigenvalue problem A * x = w * M * x

    Outputs:

    w -- An real array of k eigenvalues 

    v -- An array of k real eigenvectors, k[i] is the eigenvector corresponding
         to the eigenvector w[i]

    Optional Inputs:

    ncv -- Number of Lanczos vectors generated, ncv must be greater than k
           and is recommended to be ncv > 2*k

    which -- String specifying which eigenvectors to compute.
             Compute the k 
             'LA' - largest (algebraic) eigenvalues.
             'SA' - smallest (algebraic) eigenvalues.
             'LM' - largest (in magnitude) eigenvalues.
             'SM' - smallest (in magnitude) eigenvalues. 
             'BE' - eigenvalues, half from each end of the
                    spectrum.  When NEV is odd, compute one more from the
                    high end than from the low end.

    maxiter --  Maximum number of Arnoldi update iterations allowed

    tol -- Relative accuracy for eigenvalues (stopping criterion)

    return_eigenvectors -- True|False, return eigenvectors 

    """
    A=sb.asanyarray(A)
    n,ny=A.shape
    try:
        n==ny
    except:
        raise AttributeError("matrix is not square")
    if M is not None:
        raise NotImplementedError("generalized eigenproblem not supported yet")
    if ncv is None:
        ncv=2*k+1
    ncv=min(ncv,n)
    if maxiter==None:
        maxiter=n*10


    # guess type        
    resid = sb.zeros(n,'f')
    try:
        typ = A.dtype.char
    except AttributeError:
        typ = A.matvec(resid).dtype.char
    if typ not in 'fd':
        raise ValueError("matrix type must be 'f' or 'd'")

    # some sanity checks
    if k <= 0:
        raise ValueError("k must be positive, k=%d"%k)
    if k == n:
        raise ValueError("k must be less than rank(A), k=%d"%k)
    if maxiter <= 0:
        raise ValueError("maxiter must be positive, maxiter=%d"%maxiter)
    whiches=['LM','SM','LA','SA','BE']
    if which not in whiches:
        raise ValueError("which must be one of %s"%' '.join(whiches))
    if ncv > n or ncv < k:
        raise ValueError("ncv must be k<=ncv<=n, ncv=%s"%ncv)

    # assign solver and postprocessor
    ltr = _type_conv[typ]
    eigsolver = _arpack.__dict__[ltr+'saupd']
    eigextract = _arpack.__dict__[ltr+'seupd']
    matvec = get_matvec(A)

    v = sb.zeros((n,ncv),typ)
    resid = sb.zeros(n,typ)
    workd = sb.zeros(3*n,typ)
    workl = sb.zeros(ncv*(ncv+8),typ)
    iparam = sb.zeros(11,'int')
    ipntr = sb.zeros(11,'int')
    info = 0
    ido = 0

    # only supported mode is 1: Ax=lx
    ishfts = 1
    mode1 = 1
    bmat='I'
    iparam[0] = ishfts
    iparam[2] = maxiter
    iparam[6] = mode1


    while True:
        ido,resid,v,iparam,ipntr,info =\
            eigsolver(ido,bmat,which,k,tol,resid,v,iparam,ipntr,
               workd,workl,info)
        if (ido == -1 or ido == 1):
            xslice = slice(ipntr[0]-1, ipntr[0]-1+n)
            yslice = slice(ipntr[1]-1, ipntr[1]-1+n)
            workd[yslice]=matvec(workd[xslice])
        else:
            break

    if  info < -1 :
        raise RuntimeError("Error info=%d in arpack"%info)
        return None
    if info == -1:
        warnings.warn("Maximum number of iterations taken: %s"%iparam[2])

    # now extract eigenvalues and (optionally) eigenvectors        
    rvec = return_eigenvectors
    ierr = 0
    howmny = 'A' # return all eigenvectors
    sselect = sb.zeros(ncv,'int') # unused 
    sigma = 0.0 # no shifts, not implemented

    d,z,info =\
             eigextract(rvec,howmny,sselect,sigma,
                        bmat,which, k,tol,resid,v,iparam[0:7],ipntr,
                        workd[0:2*n],workl,ierr)

    if ierr != 0:
        raise RuntimeError("Error info=%d in arpack"%info)
        return None
    if return_eigenvectors:
        return d,z
    return d