Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/mathlibs/src/arpack/dseupd.f
5199 views
1
c\BeginDoc
2
c
3
c\Name: dseupd
4
c
5
c\Description:
6
c
7
c This subroutine returns the converged approximations to eigenvalues
8
c of A*z = lambda*B*z and (optionally):
9
c
10
c (1) the corresponding approximate eigenvectors,
11
c
12
c (2) an orthonormal (Lanczos) basis for the associated approximate
13
c invariant subspace,
14
c
15
c (3) Both.
16
c
17
c There is negligible additional cost to obtain eigenvectors. An orthonormal
18
c (Lanczos) basis is always computed. There is an additional storage cost
19
c of n*nev if both are requested (in this case a separate array Z must be
20
c supplied).
21
c
22
c These quantities are obtained from the Lanczos factorization computed
23
c by DSAUPD for the linear operator OP prescribed by the MODE selection
24
c (see IPARAM(7) in DSAUPD documentation.) DSAUPD must be called before
25
c this routine is called. These approximate eigenvalues and vectors are
26
c commonly called Ritz values and Ritz vectors respectively. They are
27
c referred to as such in the comments that follow. The computed orthonormal
28
c basis for the invariant subspace corresponding to these Ritz values is
29
c referred to as a Lanczos basis.
30
c
31
c See documentation in the header of the subroutine DSAUPD for a definition
32
c of OP as well as other terms and the relation of computed Ritz values
33
c and vectors of OP with respect to the given problem A*z = lambda*B*z.
34
c
35
c The approximate eigenvalues of the original problem are returned in
36
c ascending algebraic order. The user may elect to call this routine
37
c once for each desired Ritz vector and store it peripherally if desired.
38
c There is also the option of computing a selected set of these vectors
39
c with a single call.
40
c
41
c\Usage:
42
c call dseupd
43
c ( RVEC, HOWMNY, SELECT, D, Z, LDZ, SIGMA, BMAT, N, WHICH, NEV, TOL,
44
c RESID, NCV, V, LDV, IPARAM, IPNTR, WORKD, WORKL, LWORKL, INFO )
45
c
46
c RVEC LOGICAL (INPUT)
47
c Specifies whether Ritz vectors corresponding to the Ritz value
48
c approximations to the eigenproblem A*z = lambda*B*z are computed.
49
c
50
c RVEC = .FALSE. Compute Ritz values only.
51
c
52
c RVEC = .TRUE. Compute Ritz vectors.
53
c
54
c HOWMNY Character*1 (INPUT)
55
c Specifies how many Ritz vectors are wanted and the form of Z
56
c the matrix of Ritz vectors. See remark 1 below.
57
c = 'A': compute NEV Ritz vectors;
58
c = 'S': compute some of the Ritz vectors, specified
59
c by the logical array SELECT.
60
c
61
c SELECT Logical array of dimension NEV. (INPUT)
62
c If HOWMNY = 'S', SELECT specifies the Ritz vectors to be
63
c computed. To select the Ritz vector corresponding to a
64
c Ritz value D(j), SELECT(j) must be set to .TRUE..
65
c If HOWMNY = 'A' , SELECT is not referenced.
66
c
67
c D Double precision array of dimension NEV. (OUTPUT)
68
c On exit, D contains the Ritz value approximations to the
69
c eigenvalues of A*z = lambda*B*z. The values are returned
70
c in ascending order. If IPARAM(7) = 3,4,5 then D represents
71
c the Ritz values of OP computed by dsaupd transformed to
72
c those of the original eigensystem A*z = lambda*B*z. If
73
c IPARAM(7) = 1,2 then the Ritz values of OP are the same
74
c as the those of A*z = lambda*B*z.
75
c
76
c Z Double precision N by NEV array if HOWMNY = 'A'. (OUTPUT)
77
c On exit, Z contains the B-orthonormal Ritz vectors of the
78
c eigensystem A*z = lambda*B*z corresponding to the Ritz
79
c value approximations.
80
c If RVEC = .FALSE. then Z is not referenced.
81
c NOTE: The array Z may be set equal to first NEV columns of the
82
c Arnoldi/Lanczos basis array V computed by DSAUPD.
83
c
84
c LDZ Integer. (INPUT)
85
c The leading dimension of the array Z. If Ritz vectors are
86
c desired, then LDZ .ge. max( 1, N ). In any case, LDZ .ge. 1.
87
c
88
c SIGMA Double precision (INPUT)
89
c If IPARAM(7) = 3,4,5 represents the shift. Not referenced if
90
c IPARAM(7) = 1 or 2.
91
c
92
c
93
c **** The remaining arguments MUST be the same as for the ****
94
c **** call to DNAUPD that was just completed. ****
95
c
96
c NOTE: The remaining arguments
97
c
98
c BMAT, N, WHICH, NEV, TOL, RESID, NCV, V, LDV, IPARAM, IPNTR,
99
c WORKD, WORKL, LWORKL, INFO
100
c
101
c must be passed directly to DSEUPD following the last call
102
c to DSAUPD. These arguments MUST NOT BE MODIFIED between
103
c the the last call to DSAUPD and the call to DSEUPD.
104
c
105
c Two of these parameters (WORKL, INFO) are also output parameters:
106
c
107
c WORKL Double precision work array of length LWORKL. (OUTPUT/WORKSPACE)
108
c WORKL(1:4*ncv) contains information obtained in
109
c dsaupd. They are not changed by dseupd.
110
c WORKL(4*ncv+1:ncv*ncv+8*ncv) holds the
111
c untransformed Ritz values, the computed error estimates,
112
c and the associated eigenvector matrix of H.
113
c
114
c Note: IPNTR(8:10) contains the pointer into WORKL for addresses
115
c of the above information computed by dseupd.
116
c -------------------------------------------------------------
117
c IPNTR(8): pointer to the NCV RITZ values of the original system.
118
c IPNTR(9): pointer to the NCV corresponding error bounds.
119
c IPNTR(10): pointer to the NCV by NCV matrix of eigenvectors
120
c of the tridiagonal matrix T. Only referenced by
121
c dseupd if RVEC = .TRUE. See Remarks.
122
c -------------------------------------------------------------
123
c
124
c INFO Integer. (OUTPUT)
125
c Error flag on output.
126
c = 0: Normal exit.
127
c = -1: N must be positive.
128
c = -2: NEV must be positive.
129
c = -3: NCV must be greater than NEV and less than or equal to N.
130
c = -5: WHICH must be one of 'LM', 'SM', 'LA', 'SA' or 'BE'.
131
c = -6: BMAT must be one of 'I' or 'G'.
132
c = -7: Length of private work WORKL array is not sufficient.
133
c = -8: Error return from trid. eigenvalue calculation;
134
c Information error from LAPACK routine dsteqr.
135
c = -9: Starting vector is zero.
136
c = -10: IPARAM(7) must be 1,2,3,4,5.
137
c = -11: IPARAM(7) = 1 and BMAT = 'G' are incompatible.
138
c = -12: NEV and WHICH = 'BE' are incompatible.
139
c = -14: DSAUPD did not find any eigenvalues to sufficient
140
c accuracy.
141
c = -15: HOWMNY must be one of 'A' or 'S' if RVEC = .true.
142
c = -16: HOWMNY = 'S' not yet implemented
143
c
144
c\BeginLib
145
c
146
c\References:
147
c 1. D.C. Sorensen, "Implicit Application of Polynomial Filters in
148
c a k-Step Arnoldi Method", SIAM J. Matr. Anal. Apps., 13 (1992),
149
c pp 357-385.
150
c 2. R.B. Lehoucq, "Analysis and Implementation of an Implicitly
151
c Restarted Arnoldi Iteration", Rice University Technical Report
152
c TR95-13, Department of Computational and Applied Mathematics.
153
c 3. B.N. Parlett, "The Symmetric Eigenvalue Problem". Prentice-Hall,
154
c 1980.
155
c 4. B.N. Parlett, B. Nour-Omid, "Towards a Black Box Lanczos Program",
156
c Computer Physics Communications, 53 (1989), pp 169-179.
157
c 5. B. Nour-Omid, B.N. Parlett, T. Ericson, P.S. Jensen, "How to
158
c Implement the Spectral Transformation", Math. Comp., 48 (1987),
159
c pp 663-673.
160
c 6. R.G. Grimes, J.G. Lewis and H.D. Simon, "A Shifted Block Lanczos
161
c Algorithm for Solving Sparse Symmetric Generalized Eigenproblems",
162
c SIAM J. Matr. Anal. Apps., January (1993).
163
c 7. L. Reichel, W.B. Gragg, "Algorithm 686: FORTRAN Subroutines
164
c for Updating the QR decomposition", ACM TOMS, December 1990,
165
c Volume 16 Number 4, pp 369-377.
166
c
167
c\Remarks
168
c 1. The converged Ritz values are always returned in increasing
169
c (algebraic) order.
170
c
171
c 2. Currently only HOWMNY = 'A' is implemented. It is included at this
172
c stage for the user who wants to incorporate it.
173
c
174
c\Routines called:
175
c dsesrt ARPACK routine that sorts an array X, and applies the
176
c corresponding permutation to a matrix A.
177
c dsortr dsortr ARPACK sorting routine.
178
c ivout ARPACK utility routine that prints integers.
179
c dvout ARPACK utility routine that prints vectors.
180
c dgeqr2 LAPACK routine that computes the QR factorization of
181
c a matrix.
182
c dlacpy LAPACK matrix copy routine.
183
c dlamch LAPACK routine that determines machine constants.
184
c dorm2r LAPACK routine that applies an orthogonal matrix in
185
c factored form.
186
c dsteqr LAPACK routine that computes eigenvalues and eigenvectors
187
c of a tridiagonal matrix.
188
c dger Level 2 BLAS rank one update to a matrix.
189
c dcopy Level 1 BLAS that copies one vector to another .
190
c dnrm2 Level 1 BLAS that computes the norm of a vector.
191
c dscal Level 1 BLAS that scales a vector.
192
c dswap Level 1 BLAS that swaps the contents of two vectors.
193
194
c\Authors
195
c Danny Sorensen Phuong Vu
196
c Richard Lehoucq CRPC / Rice University
197
c Chao Yang Houston, Texas
198
c Dept. of Computational &
199
c Applied Mathematics
200
c Rice University
201
c Houston, Texas
202
c
203
c\Revision history:
204
c 12/15/93: Version ' 2.1'
205
c
206
c\SCCS Information: @(#)
207
c FILE: seupd.F SID: 2.7 DATE OF SID: 8/27/96 RELEASE: 2
208
c
209
c\EndLib
210
c
211
c-----------------------------------------------------------------------
212
subroutine dseupd (rvec, howmny, select, d, z, ldz, sigma, bmat,
213
& n, which, nev, tol, resid, ncv, v, ldv, iparam,
214
& ipntr, workd, workl, lworkl, info )
215
c
216
c %----------------------------------------------------%
217
c | Include files for debugging and timing information |
218
c %----------------------------------------------------%
219
c
220
include 'debug.h'
221
include 'stat.h'
222
c
223
c %------------------%
224
c | Scalar Arguments |
225
c %------------------%
226
c
227
character bmat, howmny, which*2
228
logical rvec, select(ncv)
229
integer info, ldz, ldv, lworkl, n, ncv, nev
230
Double precision
231
& sigma, tol
232
c
233
c %-----------------%
234
c | Array Arguments |
235
c %-----------------%
236
c
237
integer iparam(7), ipntr(11)
238
Double precision
239
& d(nev), resid(n), v(ldv,ncv), z(ldz, nev),
240
& workd(2*n), workl(lworkl)
241
c
242
c %------------%
243
c | Parameters |
244
c %------------%
245
c
246
Double precision
247
& one, zero
248
parameter (one = 1.0D+0, zero = 0.0D+0)
249
c
250
c %---------------%
251
c | Local Scalars |
252
c %---------------%
253
c
254
character type*6
255
integer bounds, ierr, ih, ihb, ihd, iq, iw, j, k,
256
& ldh, ldq, mode, msglvl, nconv, next, ritz,
257
& irz, ibd, ktrord, leftptr, rghtptr, ism, ilg
258
Double precision
259
& bnorm2, rnorm, temp, thres1, thres2, tempbnd, eps23
260
logical reord
261
c
262
c %--------------%
263
c | Local Arrays |
264
c %--------------%
265
c
266
Double precision
267
& kv(2)
268
c
269
c %----------------------%
270
c | External Subroutines |
271
c %----------------------%
272
c
273
external dcopy, dger, dgeqr2, dlacpy, dorm2r, dscal,
274
& dsesrt, dsteqr, dswap, dvout, ivout, dsortr
275
c
276
c %--------------------%
277
c | External Functions |
278
c %--------------------%
279
c
280
Double precision
281
& dnrm2, dlamch
282
external dnrm2, dlamch
283
c
284
c %---------------------%
285
c | Intrinsic Functions |
286
c %---------------------%
287
c
288
intrinsic min
289
c
290
c %-----------------------%
291
c | Executable Statements |
292
c %-----------------------%
293
c
294
c %------------------------%
295
c | Set default parameters |
296
c %------------------------%
297
c
298
msglvl = mseupd
299
mode = iparam(7)
300
nconv = iparam(5)
301
info = 0
302
c
303
c %--------------%
304
c | Quick return |
305
c %--------------%
306
c
307
if (nconv .eq. 0) go to 9000
308
ierr = 0
309
c
310
if (nconv .le. 0) ierr = -14
311
if (n .le. 0) ierr = -1
312
if (nev .le. 0) ierr = -2
313
if (ncv .le. nev .or. ncv .gt. n) ierr = -3
314
if (which .ne. 'LM' .and.
315
& which .ne. 'SM' .and.
316
& which .ne. 'LA' .and.
317
& which .ne. 'SA' .and.
318
& which .ne. 'BE') ierr = -5
319
if (bmat .ne. 'I' .and. bmat .ne. 'G') ierr = -6
320
if ( (howmny .ne. 'A' .and.
321
& howmny .ne. 'P' .and.
322
& howmny .ne. 'S') .and. rvec )
323
& ierr = -15
324
if (rvec .and. howmny .eq. 'S') ierr = -16
325
c
326
if (rvec .and. lworkl .lt. ncv**2+8*ncv) ierr = -7
327
c
328
if (mode .eq. 1 .or. mode .eq. 2) then
329
type = 'REGULR'
330
else if (mode .eq. 3 ) then
331
type = 'SHIFTI'
332
else if (mode .eq. 4 ) then
333
type = 'BUCKLE'
334
else if (mode .eq. 5 ) then
335
type = 'CAYLEY'
336
else
337
ierr = -10
338
end if
339
if (mode .eq. 1 .and. bmat .eq. 'G') ierr = -11
340
if (nev .eq. 1 .and. which .eq. 'BE') ierr = -12
341
c
342
c %------------%
343
c | Error Exit |
344
c %------------%
345
c
346
if (ierr .ne. 0) then
347
info = ierr
348
go to 9000
349
end if
350
c
351
c %-------------------------------------------------------%
352
c | Pointer into WORKL for address of H, RITZ, BOUNDS, Q |
353
c | etc... and the remaining workspace. |
354
c | Also update pointer to be used on output. |
355
c | Memory is laid out as follows: |
356
c | workl(1:2*ncv) := generated tridiagonal matrix H |
357
c | The subdiagonal is stored in workl(2:ncv). |
358
c | The dead spot is workl(1) but upon exiting |
359
c | dsaupd stores the B-norm of the last residual |
360
c | vector in workl(1). We use this !!! |
361
c | workl(2*ncv+1:2*ncv+ncv) := ritz values |
362
c | The wanted values are in the first NCONV spots. |
363
c | workl(3*ncv+1:3*ncv+ncv) := computed Ritz estimates |
364
c | The wanted values are in the first NCONV spots. |
365
c | NOTE: workl(1:4*ncv) is set by dsaupd and is not |
366
c | modified by dseupd. |
367
c %-------------------------------------------------------%
368
c
369
c %-------------------------------------------------------%
370
c | The following is used and set by dseupd. |
371
c | workl(4*ncv+1:4*ncv+ncv) := used as workspace during |
372
c | computation of the eigenvectors of H. Stores |
373
c | the diagonal of H. Upon EXIT contains the NCV |
374
c | Ritz values of the original system. The first |
375
c | NCONV spots have the wanted values. If MODE = |
376
c | 1 or 2 then will equal workl(2*ncv+1:3*ncv). |
377
c | workl(5*ncv+1:5*ncv+ncv) := used as workspace during |
378
c | computation of the eigenvectors of H. Stores |
379
c | the subdiagonal of H. Upon EXIT contains the |
380
c | NCV corresponding Ritz estimates of the |
381
c | original system. The first NCONV spots have the |
382
c | wanted values. If MODE = 1,2 then will equal |
383
c | workl(3*ncv+1:4*ncv). |
384
c | workl(6*ncv+1:6*ncv+ncv*ncv) := orthogonal Q that is |
385
c | the eigenvector matrix for H as returned by |
386
c | dsteqr. Not referenced if RVEC = .False. |
387
c | Ordering follows that of workl(4*ncv+1:5*ncv) |
388
c | workl(6*ncv+ncv*ncv+1:6*ncv+ncv*ncv+2*ncv) := |
389
c | Workspace. Needed by dsteqr and by dseupd. |
390
c | GRAND total of NCV*(NCV+8) locations. |
391
c %-------------------------------------------------------%
392
c
393
c
394
ih = ipntr(5)
395
ritz = ipntr(6)
396
bounds = ipntr(7)
397
ldh = ncv
398
ldq = ncv
399
ihd = bounds + ldh
400
ihb = ihd + ldh
401
iq = ihb + ldh
402
iw = iq + ldh*ncv
403
next = iw + 2*ncv
404
ipntr(4) = next
405
ipntr(8) = ihd
406
ipntr(9) = ihb
407
ipntr(10) = iq
408
c
409
c %----------------------------------------%
410
c | irz points to the Ritz values computed |
411
c | by _seigt before exiting _saup2. |
412
c | ibd points to the Ritz estimates |
413
c | computed by _seigt before exiting |
414
c | _saup2. |
415
c %----------------------------------------%
416
c
417
irz = ipntr(11)+ncv
418
ibd = irz+ncv
419
c
420
c
421
c %---------------------------------%
422
c | Set machine dependent constant. |
423
c %---------------------------------%
424
c
425
eps23 = dlamch('Epsilon-Machine')
426
eps23 = eps23**(2.0D+0 / 3.0D+0)
427
c
428
c %---------------------------------------%
429
c | RNORM is B-norm of the RESID(1:N). |
430
c | BNORM2 is the 2 norm of B*RESID(1:N). |
431
c | Upon exit of dsaupd WORKD(1:N) has |
432
c | B*RESID(1:N). |
433
c %---------------------------------------%
434
c
435
rnorm = workl(ih)
436
if (bmat .eq. 'I') then
437
bnorm2 = rnorm
438
else if (bmat .eq. 'G') then
439
bnorm2 = dnrm2(n, workd, 1)
440
end if
441
c
442
if (rvec) then
443
c
444
c %------------------------------------------------%
445
c | Get the converged Ritz value on the boundary. |
446
c | This value will be used to dermine whether we |
447
c | need to reorder the eigenvalues and |
448
c | eigenvectors comupted by _steqr, and is |
449
c | referred to as the "threshold" value. |
450
c | |
451
c | A Ritz value gamma is said to be a wanted |
452
c | one, if |
453
c | abs(gamma) .ge. threshold, when WHICH = 'LM'; |
454
c | abs(gamma) .le. threshold, when WHICH = 'SM'; |
455
c | gamma .ge. threshold, when WHICH = 'LA'; |
456
c | gamma .le. threshold, when WHICH = 'SA'; |
457
c | gamma .le. thres1 .or. gamma .ge. thres2 |
458
c | when WHICH = 'BE'; |
459
c | |
460
c | Note: converged Ritz values and associated |
461
c | Ritz estimates have been placed in the first |
462
c | NCONV locations in workl(ritz) and |
463
c | workl(bounds) respectively. They have been |
464
c | sorted (in _saup2) according to the WHICH |
465
c | selection criterion. (Except in the case |
466
c | WHICH = 'BE', they are sorted in an increasing |
467
c | order.) |
468
c %------------------------------------------------%
469
c
470
if (which .eq. 'LM' .or. which .eq. 'SM'
471
& .or. which .eq. 'LA' .or. which .eq. 'SA' ) then
472
c
473
thres1 = workl(ritz)
474
c
475
if (msglvl .gt. 2) then
476
call dvout(logfil, 1, [thres1], ndigit,
477
& '_seupd: Threshold eigenvalue used for re-ordering')
478
end if
479
c
480
else if (which .eq. 'BE') then
481
c
482
c %------------------------------------------------%
483
c | Ritz values returned from _saup2 have been |
484
c | sorted in increasing order. Thus two |
485
c | "threshold" values (one for the small end, one |
486
c | for the large end) are in the middle. |
487
c %------------------------------------------------%
488
c
489
ism = max(nev,nconv) / 2
490
ilg = ism + 1
491
thres1 = workl(ism)
492
thres2 = workl(ilg)
493
c
494
if (msglvl .gt. 2) then
495
kv(1) = thres1
496
kv(2) = thres2
497
call dvout(logfil, 2, kv, ndigit,
498
& '_seupd: Threshold eigenvalues used for re-ordering')
499
end if
500
c
501
end if
502
c
503
c %----------------------------------------------------------%
504
c | Check to see if all converged Ritz values appear within |
505
c | the first NCONV diagonal elements returned from _seigt. |
506
c | This is done in the following way: |
507
c | |
508
c | 1) For each Ritz value obtained from _seigt, compare it |
509
c | with the threshold Ritz value computed above to |
510
c | determine whether it is a wanted one. |
511
c | |
512
c | 2) If it is wanted, then check the corresponding Ritz |
513
c | estimate to see if it has converged. If it has, set |
514
c | correponding entry in the logical array SELECT to |
515
c | .TRUE.. |
516
c | |
517
c | If SELECT(j) = .TRUE. and j > NCONV, then there is a |
518
c | converged Ritz value that does not appear at the top of |
519
c | the diagonal matrix computed by _seigt in _saup2. |
520
c | Reordering is needed. |
521
c %----------------------------------------------------------%
522
c
523
reord = .false.
524
ktrord = 0
525
do 10 j = 0, ncv-1
526
select(j+1) = .false.
527
if (which .eq. 'LM') then
528
if (abs(workl(irz+j)) .ge. abs(thres1)) then
529
tempbnd = max( eps23, abs(workl(irz+j)) )
530
if (workl(ibd+j) .le. tol*tempbnd) then
531
select(j+1) = .true.
532
end if
533
end if
534
else if (which .eq. 'SM') then
535
if (abs(workl(irz+j)) .le. abs(thres1)) then
536
tempbnd = max( eps23, abs(workl(irz+j)) )
537
if (workl(ibd+j) .le. tol*tempbnd) then
538
select(j+1) = .true.
539
end if
540
end if
541
else if (which .eq. 'LA') then
542
if (workl(irz+j) .ge. thres1) then
543
tempbnd = max( eps23, abs(workl(irz+j)) )
544
if (workl(ibd+j) .le. tol*tempbnd) then
545
select(j+1) = .true.
546
end if
547
end if
548
else if (which .eq. 'SA') then
549
if (workl(irz+j) .le. thres1) then
550
tempbnd = max( eps23, abs(workl(irz+j)) )
551
if (workl(ibd+j) .le. tol*tempbnd) then
552
select(j+1) = .true.
553
end if
554
end if
555
else if (which .eq. 'BE') then
556
if ( workl(irz+j) .le. thres1 .or.
557
& workl(irz+j) .ge. thres2 ) then
558
tempbnd = max( eps23, abs(workl(irz+j)) )
559
if (workl(ibd+j) .le. tol*tempbnd) then
560
select(j+1) = .true.
561
end if
562
end if
563
end if
564
if (j+1 .gt. nconv ) reord = select(j+1) .or. reord
565
if (select(j+1)) ktrord = ktrord + 1
566
10 continue
567
568
c %-------------------------------------------%
569
c | If KTRORD .ne. NCONV, something is wrong. |
570
c %-------------------------------------------%
571
c
572
if (msglvl .gt. 2) then
573
call ivout(logfil, 1, [ktrord], ndigit,
574
& '_seupd: Number of specified eigenvalues')
575
call ivout(logfil, 1, [nconv], ndigit,
576
& '_seupd: Number of "converged" eigenvalues')
577
end if
578
c
579
c %-----------------------------------------------------------%
580
c | Call LAPACK routine _steqr to compute the eigenvalues and |
581
c | eigenvectors of the final symmetric tridiagonal matrix H. |
582
c | Initialize the eigenvector matrix Q to the identity. |
583
c %-----------------------------------------------------------%
584
c
585
call dcopy (ncv-1, workl(ih+1), 1, workl(ihb), 1)
586
call dcopy (ncv, workl(ih+ldh), 1, workl(ihd), 1)
587
c
588
call dsteqr ('Identity', ncv, workl(ihd), workl(ihb),
589
& workl(iq), ldq, workl(iw), ierr)
590
c
591
if (ierr .ne. 0) then
592
info = -8
593
go to 9000
594
end if
595
c
596
if (msglvl .gt. 1) then
597
call dcopy (ncv, workl(iq+ncv-1), ldq, workl(iw), 1)
598
call dvout (logfil, ncv, workl(ihd), ndigit,
599
& '_seupd: NCV Ritz values of the final H matrix')
600
call dvout (logfil, ncv, workl(iw), ndigit,
601
& '_seupd: last row of the eigenvector matrix for H')
602
end if
603
c
604
if (reord) then
605
c
606
c %---------------------------------------------%
607
c | Reordered the eigenvalues and eigenvectors |
608
c | computed by _steqr so that the "converged" |
609
c | eigenvalues appear in the first NCONV |
610
c | positions of workl(ihd), and the associated |
611
c | eigenvectors appear in the first NCONV |
612
c | columns. |
613
c %---------------------------------------------%
614
c
615
leftptr = 1
616
rghtptr = ncv
617
c
618
if (ncv .eq. 1) go to 30
619
c
620
20 if (select(leftptr)) then
621
c
622
c %-------------------------------------------%
623
c | Search, from the left, for the first Ritz |
624
c | value that has not converged. |
625
c %-------------------------------------------%
626
c
627
leftptr = leftptr + 1
628
c
629
else if ( .not. select(rghtptr)) then
630
c
631
c %----------------------------------------------%
632
c | Search, from the right, the first Ritz value |
633
c | that has converged. |
634
c %----------------------------------------------%
635
c
636
rghtptr = rghtptr - 1
637
c
638
else
639
c
640
c %----------------------------------------------%
641
c | Swap the Ritz value on the left that has not |
642
c | converged with the Ritz value on the right |
643
c | that has converged. Swap the associated |
644
c | eigenvector of the tridiagonal matrix H as |
645
c | well. |
646
c %----------------------------------------------%
647
c
648
temp = workl(ihd+leftptr-1)
649
workl(ihd+leftptr-1) = workl(ihd+rghtptr-1)
650
workl(ihd+rghtptr-1) = temp
651
call dcopy(ncv, workl(iq+ncv*(leftptr-1)), 1,
652
& workl(iw), 1)
653
call dcopy(ncv, workl(iq+ncv*(rghtptr-1)), 1,
654
& workl(iq+ncv*(leftptr-1)), 1)
655
call dcopy(ncv, workl(iw), 1,
656
& workl(iq+ncv*(rghtptr-1)), 1)
657
leftptr = leftptr + 1
658
rghtptr = rghtptr - 1
659
c
660
end if
661
c
662
if (leftptr .lt. rghtptr) go to 20
663
c
664
30 end if
665
c
666
if (msglvl .gt. 2) then
667
call dvout (logfil, ncv, workl(ihd), ndigit,
668
& '_seupd: The eigenvalues of H--reordered')
669
end if
670
c
671
c %----------------------------------------%
672
c | Load the converged Ritz values into D. |
673
c %----------------------------------------%
674
c
675
call dcopy(nconv, workl(ihd), 1, d, 1)
676
c
677
else
678
c
679
c %-----------------------------------------------------%
680
c | Ritz vectors not required. Load Ritz values into D. |
681
c %-----------------------------------------------------%
682
c
683
call dcopy (nconv, workl(ritz), 1, d, 1)
684
call dcopy (ncv, workl(ritz), 1, workl(ihd), 1)
685
c
686
end if
687
c
688
c %------------------------------------------------------------------%
689
c | Transform the Ritz values and possibly vectors and corresponding |
690
c | Ritz estimates of OP to those of A*x=lambda*B*x. The Ritz values |
691
c | (and corresponding data) are returned in ascending order. |
692
c %------------------------------------------------------------------%
693
c
694
if (type .eq. 'REGULR') then
695
c
696
c %---------------------------------------------------------%
697
c | Ascending sort of wanted Ritz values, vectors and error |
698
c | bounds. Not necessary if only Ritz values are desired. |
699
c %---------------------------------------------------------%
700
c
701
if (rvec) then
702
call dsesrt ('LA', rvec , nconv, d, ncv, workl(iq), ldq)
703
else
704
call dcopy (ncv, workl(bounds), 1, workl(ihb), 1)
705
end if
706
c
707
else
708
c
709
c %-------------------------------------------------------------%
710
c | * Make a copy of all the Ritz values. |
711
c | * Transform the Ritz values back to the original system. |
712
c | For TYPE = 'SHIFTI' the transformation is |
713
c | lambda = 1/theta + sigma |
714
c | For TYPE = 'BUCKLE' the transformation is |
715
c | lambda = sigma * theta / ( theta - 1 ) |
716
c | For TYPE = 'CAYLEY' the transformation is |
717
c | lambda = sigma * (theta + 1) / (theta - 1 ) |
718
c | where the theta are the Ritz values returned by dsaupd. |
719
c | NOTES: |
720
c | *The Ritz vectors are not affected by the transformation. |
721
c | They are only reordered. |
722
c %-------------------------------------------------------------%
723
c
724
call dcopy (ncv, workl(ihd), 1, workl(iw), 1)
725
if (type .eq. 'SHIFTI') then
726
do 40 k=1, ncv
727
workl(ihd+k-1) = one / workl(ihd+k-1) + sigma
728
40 continue
729
else if (type .eq. 'BUCKLE') then
730
do 50 k=1, ncv
731
workl(ihd+k-1) = sigma * workl(ihd+k-1) /
732
& (workl(ihd+k-1) - one)
733
50 continue
734
else if (type .eq. 'CAYLEY') then
735
do 60 k=1, ncv
736
workl(ihd+k-1) = sigma * (workl(ihd+k-1) + one) /
737
& (workl(ihd+k-1) - one)
738
60 continue
739
end if
740
c
741
c %-------------------------------------------------------------%
742
c | * Store the wanted NCONV lambda values into D. |
743
c | * Sort the NCONV wanted lambda in WORKL(IHD:IHD+NCONV-1) |
744
c | into ascending order and apply sort to the NCONV theta |
745
c | values in the transformed system. We'll need this to |
746
c | compute Ritz estimates in the original system. |
747
c | * Finally sort the lambda's into ascending order and apply |
748
c | to Ritz vectors if wanted. Else just sort lambda's into |
749
c | ascending order. |
750
c | NOTES: |
751
c | *workl(iw:iw+ncv-1) contain the theta ordered so that they |
752
c | match the ordering of the lambda. We'll use them again for |
753
c | Ritz vector purification. |
754
c %-------------------------------------------------------------%
755
c
756
call dcopy (nconv, workl(ihd), 1, d, 1)
757
call dsortr ('LA', .true., nconv, workl(ihd), workl(iw))
758
if (rvec) then
759
call dsesrt ('LA', rvec , nconv, d, ncv, workl(iq), ldq)
760
else
761
call dcopy (ncv, workl(bounds), 1, workl(ihb), 1)
762
call dscal (ncv, bnorm2/rnorm, workl(ihb), 1)
763
call dsortr ('LA', .true., nconv, d, workl(ihb))
764
end if
765
c
766
end if
767
c
768
c %------------------------------------------------%
769
c | Compute the Ritz vectors. Transform the wanted |
770
c | eigenvectors of the symmetric tridiagonal H by |
771
c | the Lanczos basis matrix V. |
772
c %------------------------------------------------%
773
c
774
if (rvec .and. howmny .eq. 'A') then
775
c
776
c %----------------------------------------------------------%
777
c | Compute the QR factorization of the matrix representing |
778
c | the wanted invariant subspace located in the first NCONV |
779
c | columns of workl(iq,ldq). |
780
c %----------------------------------------------------------%
781
c
782
call dgeqr2 (ncv, nconv, workl(iq), ldq, workl(iw+ncv),
783
& workl(ihb), ierr)
784
c
785
c
786
c %--------------------------------------------------------%
787
c | * Postmultiply V by Q. |
788
c | * Copy the first NCONV columns of VQ into Z. |
789
c | The N by NCONV matrix Z is now a matrix representation |
790
c | of the approximate invariant subspace associated with |
791
c | the Ritz values in workl(ihd). |
792
c %--------------------------------------------------------%
793
c
794
call dorm2r ('Right', 'Notranspose', n, ncv, nconv, workl(iq),
795
& ldq, workl(iw+ncv), v, ldv, workd(n+1), ierr)
796
call dlacpy ('All', n, nconv, v, ldv, z, ldz)
797
c
798
c %-----------------------------------------------------%
799
c | In order to compute the Ritz estimates for the Ritz |
800
c | values in both systems, need the last row of the |
801
c | eigenvector matrix. Remember, it's in factored form |
802
c %-----------------------------------------------------%
803
c
804
do 65 j = 1, ncv-1
805
workl(ihb+j-1) = zero
806
65 continue
807
workl(ihb+ncv-1) = one
808
call dorm2r ('Left', 'Transpose', ncv, 1, nconv, workl(iq),
809
& ldq, workl(iw+ncv), workl(ihb), ncv, temp, ierr)
810
c
811
else if (rvec .and. howmny .eq. 'S') then
812
c
813
c Not yet implemented. See remark 2 above.
814
c
815
end if
816
c
817
if (type .eq. 'REGULR' .and. rvec) then
818
c
819
do 70 j=1, ncv
820
workl(ihb+j-1) = rnorm * abs( workl(ihb+j-1) )
821
70 continue
822
c
823
else if (type .ne. 'REGULR' .and. rvec) then
824
c
825
c %-------------------------------------------------%
826
c | * Determine Ritz estimates of the theta. |
827
c | If RVEC = .true. then compute Ritz estimates |
828
c | of the theta. |
829
c | If RVEC = .false. then copy Ritz estimates |
830
c | as computed by dsaupd. |
831
c | * Determine Ritz estimates of the lambda. |
832
c %-------------------------------------------------%
833
c
834
call dscal (ncv, bnorm2, workl(ihb), 1)
835
if (type .eq. 'SHIFTI') then
836
c
837
do 80 k=1, ncv
838
workl(ihb+k-1) = abs( workl(ihb+k-1) ) / workl(iw+k-1)**2
839
80 continue
840
c
841
else if (type .eq. 'BUCKLE') then
842
c
843
do 90 k=1, ncv
844
workl(ihb+k-1) = sigma * abs( workl(ihb+k-1) ) /
845
& ( workl(iw+k-1)-one )**2
846
90 continue
847
c
848
else if (type .eq. 'CAYLEY') then
849
c
850
do 100 k=1, ncv
851
workl(ihb+k-1) = abs( workl(ihb+k-1) /
852
& workl(iw+k-1)*(workl(iw+k-1)-one) )
853
100 continue
854
c
855
end if
856
c
857
end if
858
c
859
if (type .ne. 'REGULR' .and. msglvl .gt. 1) then
860
call dvout (logfil, nconv, d, ndigit,
861
& '_seupd: Untransformed converged Ritz values')
862
call dvout (logfil, nconv, workl(ihb), ndigit,
863
& '_seupd: Ritz estimates of the untransformed Ritz values')
864
else if (msglvl .gt. 1) then
865
call dvout (logfil, nconv, d, ndigit,
866
& '_seupd: Converged Ritz values')
867
call dvout (logfil, nconv, workl(ihb), ndigit,
868
& '_seupd: Associated Ritz estimates')
869
end if
870
c
871
c %-------------------------------------------------%
872
c | Ritz vector purification step. Formally perform |
873
c | one of inverse subspace iteration. Only used |
874
c | for MODE = 3,4,5. See reference 7 |
875
c %-------------------------------------------------%
876
c
877
if (rvec .and. (type .eq. 'SHIFTI' .or. type .eq. 'CAYLEY')) then
878
c
879
do 110 k=0, nconv-1
880
workl(iw+k) = workl(iq+k*ldq+ncv-1) / workl(iw+k)
881
110 continue
882
c
883
else if (rvec .and. type .eq. 'BUCKLE') then
884
c
885
do 120 k=0, nconv-1
886
workl(iw+k) = workl(iq+k*ldq+ncv-1) / (workl(iw+k)-one)
887
120 continue
888
c
889
end if
890
c
891
if (type .ne. 'REGULR')
892
& call dger (n, nconv, one, resid, 1, workl(iw), 1, z, ldz)
893
c
894
9000 continue
895
c
896
return
897
c
898
c %---------------%
899
c | End of dseupd |
900
c %---------------%
901
c
902
end
903
904