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