Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/mathlibs/src/arpack/dneupd.f
5213 views
1
c\BeginDoc
2
c
3
c\Name: dneupd
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 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 basis is always computed. There is an additional storage cost of n*nev
19
c if both are requested (in this case a separate array Z must be supplied).
20
c
21
c The approximate eigenvalues and eigenvectors of A*z = lambda*B*z
22
c are derived from approximate eigenvalues and eigenvectors of
23
c of the linear operator OP prescribed by the MODE selection in the
24
c call to DNAUPD. DNAUPD must be called before this routine is called.
25
c These approximate eigenvalues and vectors are commonly called Ritz
26
c values and Ritz vectors respectively. They are referred to as such
27
c in the comments that follow. The computed orthonormal basis for the
28
c invariant subspace corresponding to these Ritz values is referred to as a
29
c Schur basis.
30
c
31
c See documentation in the header of the subroutine DNAUPD for
32
c definition of OP as well as other terms and the relation of computed
33
c Ritz values and Ritz vectors of OP with respect to the given problem
34
c A*z = lambda*B*z. For a brief description, see definitions of
35
c IPARAM(7), MODE and WHICH in the documentation of DNAUPD.
36
c
37
c\Usage:
38
c call dneupd
39
c ( RVEC, HOWMNY, SELECT, DR, DI, Z, LDZ, SIGMAR, SIGMAI, WORKEV, BMAT,
40
c N, WHICH, NEV, TOL, RESID, NCV, V, LDV, IPARAM, IPNTR, WORKD, WORKL,
41
c LWORKL, INFO )
42
c
43
c\Arguments:
44
c RVEC LOGICAL (INPUT)
45
c Specifies whether a basis for the invariant subspace corresponding
46
c to the converged Ritz value approximations for the eigenproblem
47
c A*z = lambda*B*z is computed.
48
c
49
c RVEC = .FALSE. Compute Ritz values only.
50
c
51
c RVEC = .TRUE. Compute the Ritz vectors or Schur vectors.
52
c See Remarks below.
53
c
54
c HOWMNY Character*1 (INPUT)
55
c Specifies the form of the basis for the invariant subspace
56
c corresponding to the converged Ritz values that is to be computed.
57
c
58
c = 'A': Compute NEV Ritz vectors;
59
c = 'P': Compute NEV Schur vectors;
60
c = 'S': compute some of the Ritz vectors, specified
61
c by the logical array SELECT.
62
c
63
c SELECT Logical array of dimension NCV. (INPUT)
64
c If HOWMNY = 'S', SELECT specifies the Ritz vectors to be
65
c computed. To select the Ritz vector corresponding to a
66
c Ritz value (DR(j), DI(j)), SELECT(j) must be set to .TRUE..
67
c If HOWMNY = 'A' or 'P', SELECT is used as internal workspace.
68
c
69
c DR Double precision array of dimension NEV+1. (OUTPUT)
70
c If IPARAM(7) = 1,2 or 3 and SIGMAI=0.0 then on exit: DR contains
71
c the real part of the Ritz approximations to the eigenvalues of
72
c A*z = lambda*B*z.
73
c If IPARAM(7) = 3, 4 and SIGMAI is not equal to zero, then on exit:
74
c DR contains the real part of the Ritz values of OP computed by
75
c DNAUPD. A further computation must be performed by the user
76
c to transform the Ritz values computed for OP by DNAUPD to those
77
c of the original system A*z = lambda*B*z. See remark 3 below.
78
c
79
c DI Double precision array of dimension NEV+1. (OUTPUT)
80
c On exit, DI contains the imaginary part of the Ritz value
81
c approximations to the eigenvalues of A*z = lambda*B*z associated
82
c with DR.
83
c
84
c NOTE: When Ritz values are complex, they will come in complex
85
c conjugate pairs. If eigenvectors are requested, the
86
c corresponding Ritz vectors will also come in conjugate
87
c pairs and the real and imaginary parts of these are
88
c represented in two consecutive columns of the array Z
89
c (see below).
90
c
91
c Z Double precision N by NEV+1 array if RVEC = .TRUE. and HOWMNY = 'A'. (OUTPUT)
92
c On exit, if RVEC = .TRUE. and HOWMNY = 'A', then the columns of
93
c Z represent approximate eigenvectors (Ritz vectors) corresponding
94
c to the NCONV=IPARAM(5) Ritz values for eigensystem
95
c A*z = lambda*B*z.
96
c
97
c The complex Ritz vector associated with the Ritz value
98
c with positive imaginary part is stored in two consecutive
99
c columns. The first column holds the real part of the Ritz
100
c vector and the second column holds the imaginary part. The
101
c Ritz vector associated with the Ritz value with negative
102
c imaginary part is simply the complex conjugate of the Ritz vector
103
c associated with the positive imaginary part.
104
c
105
c If RVEC = .FALSE. or HOWMNY = 'P', then Z is not referenced.
106
c
107
c NOTE: If if RVEC = .TRUE. and a Schur basis is not required,
108
c the array Z may be set equal to first NEV+1 columns of the Arnoldi
109
c basis array V computed by DNAUPD. In this case the Arnoldi basis
110
c will be destroyed and overwritten with the eigenvector basis.
111
c
112
c LDZ Integer. (INPUT)
113
c The leading dimension of the array Z. If Ritz vectors are
114
c desired, then LDZ >= max( 1, N ). In any case, LDZ >= 1.
115
c
116
c SIGMAR Double precision (INPUT)
117
c If IPARAM(7) = 3 or 4, represents the real part of the shift.
118
c Not referenced if IPARAM(7) = 1 or 2.
119
c
120
c SIGMAI Double precision (INPUT)
121
c If IPARAM(7) = 3 or 4, represents the imaginary part of the shift.
122
c Not referenced if IPARAM(7) = 1 or 2. See remark 3 below.
123
c
124
c WORKEV Double precision work array of dimension 3*NCV. (WORKSPACE)
125
c
126
c **** The remaining arguments MUST be the same as for the ****
127
c **** call to DNAUPD that was just completed. ****
128
c
129
c NOTE: The remaining arguments
130
c
131
c BMAT, N, WHICH, NEV, TOL, RESID, NCV, V, LDV, IPARAM, IPNTR,
132
c WORKD, WORKL, LWORKL, INFO
133
c
134
c must be passed directly to DNEUPD following the last call
135
c to DNAUPD. These arguments MUST NOT BE MODIFIED between
136
c the the last call to DNAUPD and the call to DNEUPD.
137
c
138
c Three of these parameters (V, WORKL, INFO) are also output parameters:
139
c
140
c V Double precision N by NCV array. (INPUT/OUTPUT)
141
c
142
c Upon INPUT: the NCV columns of V contain the Arnoldi basis
143
c vectors for OP as constructed by DNAUPD .
144
c
145
c Upon OUTPUT: If RVEC = .TRUE. the first NCONV=IPARAM(5) columns
146
c contain approximate Schur vectors that span the
147
c desired invariant subspace. See Remark 2 below.
148
c
149
c NOTE: If the array Z has been set equal to first NEV+1 columns
150
c of the array V and RVEC=.TRUE. and HOWMNY= 'A', then the
151
c Arnoldi basis held by V has been overwritten by the desired
152
c Ritz vectors. If a separate array Z has been passed then
153
c the first NCONV=IPARAM(5) columns of V will contain approximate
154
c Schur vectors that span the desired invariant subspace.
155
c
156
c WORKL Double precision work array of length LWORKL. (OUTPUT/WORKSPACE)
157
c WORKL(1:ncv*ncv+3*ncv) contains information obtained in
158
c dnaupd. They are not changed by dneupd.
159
c WORKL(ncv*ncv+3*ncv+1:3*ncv*ncv+6*ncv) holds the
160
c real and imaginary part of the untransformed Ritz values,
161
c the upper quasi-triangular matrix for H, and the
162
c associated matrix representation of the invariant subspace for H.
163
c
164
c Note: IPNTR(9:13) contains the pointer into WORKL for addresses
165
c of the above information computed by dneupd.
166
c -------------------------------------------------------------
167
c IPNTR(9): pointer to the real part of the NCV RITZ values of the
168
c original system.
169
c IPNTR(10): pointer to the imaginary part of the NCV RITZ values of
170
c the original system.
171
c IPNTR(11): pointer to the NCV corresponding error bounds.
172
c IPNTR(12): pointer to the NCV by NCV upper quasi-triangular
173
c Schur matrix for H.
174
c IPNTR(13): pointer to the NCV by NCV matrix of eigenvectors
175
c of the upper Hessenberg matrix H. Only referenced by
176
c dneupd if RVEC = .TRUE. See Remark 2 below.
177
c -------------------------------------------------------------
178
c
179
c INFO Integer. (OUTPUT)
180
c Error flag on output.
181
c
182
c = 0: Normal exit.
183
c
184
c = 1: The Schur form computed by LAPACK routine dlahqr
185
c could not be reordered by LAPACK routine dtrsen.
186
c Re-enter subroutine dneupd with IPARAM(5)=NCV and
187
c increase the size of the arrays DR and DI to have
188
c dimension at least dimension NCV and allocate at least NCV
189
c columns for Z. NOTE: Not necessary if Z and V share
190
c the same space. Please notify the authors if this error
191
c occurs.
192
c
193
c = -1: N must be positive.
194
c = -2: NEV must be positive.
195
c = -3: NCV-NEV >= 2 and less than or equal to N.
196
c = -5: WHICH must be one of 'LM', 'SM', 'LR', 'SR', 'LI', 'SI'
197
c = -6: BMAT must be one of 'I' or 'G'.
198
c = -7: Length of private work WORKL array is not sufficient.
199
c = -8: Error return from calculation of a real Schur form.
200
c Informational error from LAPACK routine dlahqr.
201
c = -9: Error return from calculation of eigenvectors.
202
c Informational error from LAPACK routine dtrevc.
203
c = -10: IPARAM(7) must be 1,2,3,4.
204
c = -11: IPARAM(7) = 1 and BMAT = 'G' are incompatible.
205
c = -12: HOWMNY = 'S' not yet implemented
206
c = -13: HOWMNY must be one of 'A' or 'P' if RVEC = .true.
207
c = -14: DNAUPD did not find any eigenvalues to sufficient
208
c accuracy.
209
c
210
c\BeginLib
211
c
212
c\References:
213
c 1. D.C. Sorensen, "Implicit Application of Polynomial Filters in
214
c a k-Step Arnoldi Method", SIAM J. Matr. Anal. Apps., 13 (1992),
215
c pp 357-385.
216
c 2. R.B. Lehoucq, "Analysis and Implementation of an Implicitly
217
c Restarted Arnoldi Iteration", Rice University Technical Report
218
c TR95-13, Department of Computational and Applied Mathematics.
219
c 3. B.N. Parlett & Y. Saad, "Complex Shift and Invert Strategies for
220
c Real Matrices", Linear Algebra and its Applications, vol 88/89,
221
c pp 575-595, (1987).
222
c
223
c\Routines called:
224
c ivout ARPACK utility routine that prints integers.
225
c dmout ARPACK utility routine that prints matrices
226
c dvout ARPACK utility routine that prints vectors.
227
c dgeqr2 LAPACK routine that computes the QR factorization of
228
c a matrix.
229
c dlacpy LAPACK matrix copy routine.
230
c dlahqr LAPACK routine to compute the real Schur form of an
231
c upper Hessenberg matrix.
232
c dlamch LAPACK routine that determines machine constants.
233
c dlapy2 LAPACK routine to compute sqrt(x**2+y**2) carefully.
234
c dlaset LAPACK matrix initialization routine.
235
c dorm2r LAPACK routine that applies an orthogonal matrix in
236
c factored form.
237
c dtrevc LAPACK routine to compute the eigenvectors of a matrix
238
c in upper quasi-triangular form.
239
c dtrsen LAPACK routine that re-orders the Schur form.
240
c dtrmm Level 3 BLAS matrix times an upper triangular matrix.
241
c dger Level 2 BLAS rank one update to a matrix.
242
c dcopy Level 1 BLAS that copies one vector to another .
243
c ddot Level 1 BLAS that computes the scalar product of two vectors.
244
c dnrm2 Level 1 BLAS that computes the norm of a vector.
245
c dscal Level 1 BLAS that scales a vector.
246
c
247
c\Remarks
248
c
249
c 1. Currently only HOWMNY = 'A' and 'P' are implemented.
250
c
251
c Let X' denote the transpose of X.
252
c
253
c 2. Schur vectors are an orthogonal representation for the basis of
254
c Ritz vectors. Thus, their numerical properties are often superior.
255
c If RVEC = .TRUE. then the relationship
256
c A * V(:,1:IPARAM(5)) = V(:,1:IPARAM(5)) * T, and
257
c V(:,1:IPARAM(5))' * V(:,1:IPARAM(5)) = I are approximately satisfied.
258
c Here T is the leading submatrix of order IPARAM(5) of the real
259
c upper quasi-triangular matrix stored workl(ipntr(12)). That is,
260
c T is block upper triangular with 1-by-1 and 2-by-2 diagonal blocks;
261
c each 2-by-2 diagonal block has its diagonal elements equal and its
262
c off-diagonal elements of opposite sign. Corresponding to each 2-by-2
263
c diagonal block is a complex conjugate pair of Ritz values. The real
264
c Ritz values are stored on the diagonal of T.
265
c
266
c 3. If IPARAM(7) = 3 or 4 and SIGMAI is not equal zero, then the user must
267
c form the IPARAM(5) Rayleigh quotients in order to transform the Ritz
268
c values computed by DNAUPD for OP to those of A*z = lambda*B*z.
269
c Set RVEC = .true. and HOWMNY = 'A', and
270
c compute
271
c Z(:,I)' * A * Z(:,I) if DI(I) = 0.
272
c If DI(I) is not equal to zero and DI(I+1) = - D(I),
273
c then the desired real and imaginary parts of the Ritz value are
274
c Z(:,I)' * A * Z(:,I) + Z(:,I+1)' * A * Z(:,I+1),
275
c Z(:,I)' * A * Z(:,I+1) - Z(:,I+1)' * A * Z(:,I), respectively.
276
c Another possibility is to set RVEC = .true. and HOWMNY = 'P' and
277
c compute V(:,1:IPARAM(5))' * A * V(:,1:IPARAM(5)) and then an upper
278
c quasi-triangular matrix of order IPARAM(5) is computed. See remark
279
c 2 above.
280
c
281
c\Authors
282
c Danny Sorensen Phuong Vu
283
c Richard Lehoucq CRPC / Rice University
284
c Chao Yang Houston, Texas
285
c Dept. of Computational &
286
c Applied Mathematics
287
c Rice University
288
c Houston, Texas
289
c
290
c\SCCS Information: @(#)
291
c FILE: neupd.F SID: 2.5 DATE OF SID: 7/31/96 RELEASE: 2
292
c
293
c\EndLib
294
c
295
c-----------------------------------------------------------------------
296
subroutine dneupd (rvec, howmny, select, dr, di, z, ldz, sigmar,
297
& sigmai, workev, bmat, n, which, nev, tol,
298
& resid, ncv, v, ldv, iparam, ipntr, workd,
299
& workl, lworkl, info)
300
c
301
c %----------------------------------------------------%
302
c | Include files for debugging and timing information |
303
c %----------------------------------------------------%
304
c
305
include 'debug.h'
306
include 'stat.h'
307
c
308
c %------------------%
309
c | Scalar Arguments |
310
c %------------------%
311
c
312
character bmat, howmny, which*2
313
logical rvec
314
integer info, ldz, ldv, lworkl, n, ncv, nev
315
Double precision
316
& sigmar, sigmai, tol
317
c
318
c %-----------------%
319
c | Array Arguments |
320
c %-----------------%
321
c
322
integer iparam(11), ipntr(14)
323
logical select(ncv)
324
Double precision
325
& dr(nev+1), di(nev+1), resid(n), v(ldv,ncv), z(ldz,*),
326
& workd(3*n), workl(lworkl), workev(3*ncv)
327
c
328
c %------------%
329
c | Parameters |
330
c %------------%
331
c
332
Double precision
333
& one, zero
334
parameter (one = 1.0D+0, zero = 0.0D+0)
335
c
336
c %---------------%
337
c | Local Scalars |
338
c %---------------%
339
c
340
character type*6
341
integer bounds, ierr, ih, ihbds, iheigr, iheigi, iconj, nconv,
342
& invsub, iuptri, iwev, iwork(1), j, k, ktrord,
343
& ldh, ldq, mode, msglvl, outncv, ritzr, ritzi, wri, wrr,
344
& irr, iri, ibd
345
logical reord
346
Double precision
347
& conds, rnorm, sep, temp, thres, vl(1,1), temp1, eps23
348
c
349
c %----------------------%
350
c | External Subroutines |
351
c %----------------------%
352
c
353
external dcopy, dger, dgeqr2, dlacpy, dlahqr, dlaset, dmout,
354
& dorm2r, dtrevc, dtrmm, dtrsen, dscal, dvout, ivout
355
c
356
c %--------------------%
357
c | External Functions |
358
c %--------------------%
359
c
360
Double precision
361
& dlapy2, dnrm2, dlamch, ddot
362
external dlapy2, dnrm2, dlamch, ddot
363
c
364
c %---------------------%
365
c | Intrinsic Functions |
366
c %---------------------%
367
c
368
intrinsic abs, min, sqrt
369
c
370
c %-----------------------%
371
c | Executable Statements |
372
c %-----------------------%
373
c
374
c %------------------------%
375
c | Set default parameters |
376
c %------------------------%
377
c
378
msglvl = mneupd
379
mode = iparam(7)
380
nconv = iparam(5)
381
info = 0
382
c
383
c %---------------------------------%
384
c | Get machine dependent constant. |
385
c %---------------------------------%
386
c
387
eps23 = dlamch('Epsilon-Machine')
388
eps23 = eps23**(2.0D+0 / 3.0D+0)
389
c
390
c %--------------%
391
c | Quick return |
392
c %--------------%
393
c
394
ierr = 0
395
c
396
if (nconv .le. 0) then
397
ierr = -14
398
else if (n .le. 0) then
399
ierr = -1
400
else if (nev .le. 0) then
401
ierr = -2
402
else if (ncv .le. nev+1 .or. ncv .gt. n) then
403
ierr = -3
404
else if (which .ne. 'LM' .and.
405
& which .ne. 'SM' .and.
406
& which .ne. 'LR' .and.
407
& which .ne. 'SR' .and.
408
& which .ne. 'LI' .and.
409
& which .ne. 'SI') then
410
ierr = -5
411
else if (bmat .ne. 'I' .and. bmat .ne. 'G') then
412
ierr = -6
413
else if (lworkl .lt. 3*ncv**2 + 6*ncv) then
414
ierr = -7
415
else if ( (howmny .ne. 'A' .and.
416
& howmny .ne. 'P' .and.
417
& howmny .ne. 'S') .and. rvec ) then
418
ierr = -13
419
else if (howmny .eq. 'S' ) then
420
ierr = -12
421
end if
422
c
423
if (mode .eq. 1 .or. mode .eq. 2) then
424
type = 'REGULR'
425
else if (mode .eq. 3 .and. sigmai .eq. zero) then
426
type = 'SHIFTI'
427
else if (mode .eq. 3 ) then
428
type = 'REALPT'
429
else if (mode .eq. 4 ) then
430
type = 'IMAGPT'
431
else
432
ierr = -10
433
end if
434
if (mode .eq. 1 .and. bmat .eq. 'G') ierr = -11
435
c
436
c %------------%
437
c | Error Exit |
438
c %------------%
439
c
440
if (ierr .ne. 0) then
441
info = ierr
442
go to 9000
443
end if
444
c
445
c %--------------------------------------------------------%
446
c | Pointer into WORKL for address of H, RITZ, BOUNDS, Q |
447
c | etc... and the remaining workspace. |
448
c | Also update pointer to be used on output. |
449
c | Memory is laid out as follows: |
450
c | workl(1:ncv*ncv) := generated Hessenberg matrix |
451
c | workl(ncv*ncv+1:ncv*ncv+2*ncv) := real and imaginary |
452
c | parts of ritz values |
453
c | workl(ncv*ncv+2*ncv+1:ncv*ncv+3*ncv) := error bounds |
454
c %--------------------------------------------------------%
455
c
456
c %-----------------------------------------------------------%
457
c | The following is used and set by DNEUPD. |
458
c | workl(ncv*ncv+3*ncv+1:ncv*ncv+4*ncv) := The untransformed |
459
c | real part of the Ritz values. |
460
c | workl(ncv*ncv+4*ncv+1:ncv*ncv+5*ncv) := The untransformed |
461
c | imaginary part of the Ritz values. |
462
c | workl(ncv*ncv+5*ncv+1:ncv*ncv+6*ncv) := The untransformed |
463
c | error bounds of the Ritz values |
464
c | workl(ncv*ncv+6*ncv+1:2*ncv*ncv+6*ncv) := Holds the upper |
465
c | quasi-triangular matrix for H |
466
c | workl(2*ncv*ncv+6*ncv+1: 3*ncv*ncv+6*ncv) := Holds the |
467
c | associated matrix representation of the invariant |
468
c | subspace for H. |
469
c | GRAND total of NCV * ( 3 * NCV + 6 ) locations. |
470
c %-----------------------------------------------------------%
471
c
472
ih = ipntr(5)
473
ritzr = ipntr(6)
474
ritzi = ipntr(7)
475
bounds = ipntr(8)
476
ldh = ncv
477
ldq = ncv
478
iheigr = bounds + ldh
479
iheigi = iheigr + ldh
480
ihbds = iheigi + ldh
481
iuptri = ihbds + ldh
482
invsub = iuptri + ldh*ncv
483
ipntr(9) = iheigr
484
ipntr(10) = iheigi
485
ipntr(11) = ihbds
486
ipntr(12) = iuptri
487
ipntr(13) = invsub
488
wrr = 1
489
wri = ncv + 1
490
iwev = wri + ncv
491
c
492
c %-----------------------------------------%
493
c | irr points to the REAL part of the Ritz |
494
c | values computed by _neigh before |
495
c | exiting _naup2. |
496
c | iri points to the IMAGINARY part of the |
497
c | Ritz values computed by _neigh |
498
c | before exiting _naup2. |
499
c | ibd points to the Ritz estimates |
500
c | computed by _neigh before exiting |
501
c | _naup2. |
502
c %-----------------------------------------%
503
c
504
irr = ipntr(14)+ncv*ncv
505
iri = irr+ncv
506
ibd = iri+ncv
507
c
508
c %------------------------------------%
509
c | RNORM is B-norm of the RESID(1:N). |
510
c %------------------------------------%
511
c
512
rnorm = workl(ih+2)
513
workl(ih+2) = zero
514
c
515
if (rvec) then
516
c
517
c %-------------------------------------------%
518
c | Get converged Ritz value on the boundary. |
519
c | Note: converged Ritz values have been |
520
c | placed in the first NCONV locations in |
521
c | workl(ritzr) and workl(ritzi). They have |
522
c | been sorted (in _naup2) according to the |
523
c | WHICH selection criterion. |
524
c %-------------------------------------------%
525
c
526
if (which .eq. 'LM' .or. which .eq. 'SM') then
527
thres = dlapy2( workl(ritzr), workl(ritzi) )
528
else if (which .eq. 'LR' .or. which .eq. 'SR') then
529
thres = workl(ritzr)
530
else if (which .eq. 'LI' .or. which .eq. 'SI') then
531
thres = abs( workl(ritzi) )
532
end if
533
c
534
if (msglvl .gt. 2) then
535
call dvout(logfil, 1, [thres], ndigit,
536
& '_neupd: Threshold eigenvalue used for re-ordering')
537
end if
538
c
539
c %----------------------------------------------------------%
540
c | Check to see if all converged Ritz values appear at the |
541
c | top of the upper quasi-triangular matrix computed by |
542
c | _neigh in _naup2. This is done in the following way: |
543
c | |
544
c | 1) For each Ritz value obtained from _neigh, compare it |
545
c | with the threshold Ritz value computed above to |
546
c | determine whether it is a wanted one. |
547
c | |
548
c | 2) If it is wanted, then check the corresponding Ritz |
549
c | estimate to see if it has converged. If it has, set |
550
c | correponding entry in the logical array SELECT to |
551
c | .TRUE.. |
552
c | |
553
c | If SELECT(j) = .TRUE. and j > NCONV, then there is a |
554
c | converged Ritz value that does not appear at the top of |
555
c | the upper quasi-triangular matrix computed by _neigh in |
556
c | _naup2. Reordering is needed. |
557
c %----------------------------------------------------------%
558
c
559
reord = .false.
560
ktrord = 0
561
do 10 j = 0, ncv-1
562
select(j+1) = .false.
563
if (which .eq. 'LM') then
564
if (dlapy2(workl(irr+j), workl(iri+j))
565
& .ge. thres) then
566
temp1 = max( eps23,
567
& dlapy2( workl(irr+j), workl(iri+j) ) )
568
if (workl(ibd+j) .le. tol*temp1)
569
& select(j+1) = .true.
570
end if
571
else if (which .eq. 'SM') then
572
if (dlapy2(workl(irr+j), workl(iri+j))
573
& .le. thres) then
574
temp1 = max( eps23,
575
& dlapy2( workl(irr+j), workl(iri+j) ) )
576
if (workl(ibd+j) .le. tol*temp1)
577
& select(j+1) = .true.
578
end if
579
else if (which .eq. 'LR') then
580
if (workl(irr+j) .ge. thres) then
581
temp1 = max( eps23,
582
& dlapy2( workl(irr+j), workl(iri+j) ) )
583
if (workl(ibd+j) .le. tol*temp1)
584
& select(j+1) = .true.
585
end if
586
else if (which .eq. 'SR') then
587
if (workl(irr+j) .le. thres) then
588
temp1 = max( eps23,
589
& dlapy2( workl(irr+j), workl(iri+j) ) )
590
if (workl(ibd+j) .le. tol*temp1)
591
& select(j+1) = .true.
592
end if
593
else if (which .eq. 'LI') then
594
if (abs(workl(iri+j)) .ge. thres) then
595
temp1 = max( eps23,
596
& dlapy2( workl(irr+j), workl(iri+j) ) )
597
if (workl(ibd+j) .le. tol*temp1)
598
& select(j+1) = .true.
599
end if
600
else if (which .eq. 'SI') then
601
if (abs(workl(iri+j)) .le. thres) then
602
temp1 = max( eps23,
603
& dlapy2( workl(irr+j), workl(iri+j) ) )
604
if (workl(ibd+j) .le. tol*temp1)
605
& select(j+1) = .true.
606
end if
607
end if
608
if (j+1 .gt. nconv ) reord = ( select(j+1) .or. reord )
609
if (select(j+1)) ktrord = ktrord + 1
610
10 continue
611
c
612
if (msglvl .gt. 2) then
613
call ivout(logfil, 1, [ktrord], ndigit,
614
& '_neupd: Number of specified eigenvalues')
615
call ivout(logfil, 1, [nconv], ndigit,
616
& '_neupd: Number of "converged" eigenvalues')
617
end if
618
c
619
c %-----------------------------------------------------------%
620
c | Call LAPACK routine dlahqr to compute the real Schur form |
621
c | of the upper Hessenberg matrix returned by DNAUPD. |
622
c | Make a copy of the upper Hessenberg matrix. |
623
c | Initialize the Schur vector matrix Q to the identity. |
624
c %-----------------------------------------------------------%
625
c
626
call dcopy (ldh*ncv, workl(ih), 1, workl(iuptri), 1)
627
call dlaset ('All', ncv, ncv, zero, one, workl(invsub), ldq)
628
call dlahqr (.true., .true., ncv, 1, ncv, workl(iuptri), ldh,
629
& workl(iheigr), workl(iheigi), 1, ncv,
630
& workl(invsub), ldq, ierr)
631
call dcopy (ncv, workl(invsub+ncv-1), ldq, workl(ihbds), 1)
632
c
633
if (ierr .ne. 0) then
634
info = -8
635
go to 9000
636
end if
637
c
638
if (msglvl .gt. 1) then
639
call dvout (logfil, ncv, workl(iheigr), ndigit,
640
& '_neupd: Real part of the eigenvalues of H')
641
call dvout (logfil, ncv, workl(iheigi), ndigit,
642
& '_neupd: Imaginary part of the Eigenvalues of H')
643
call dvout (logfil, ncv, workl(ihbds), ndigit,
644
& '_neupd: Last row of the Schur vector matrix')
645
if (msglvl .gt. 3) then
646
call dmout (logfil, ncv, ncv, workl(iuptri), ldh, ndigit,
647
& '_neupd: The upper quasi-triangular matrix ')
648
end if
649
end if
650
c
651
if (reord) then
652
c
653
c %-----------------------------------------------------%
654
c | Reorder the computed upper quasi-triangular matrix. |
655
c %-----------------------------------------------------%
656
c
657
call dtrsen ('None', 'V', select, ncv, workl(iuptri), ldh,
658
& workl(invsub), ldq, workl(iheigr), workl(iheigi),
659
& nconv, conds, sep, workl(ihbds), ncv, iwork, 1, ierr)
660
c
661
if (ierr .eq. 1) then
662
info = 1
663
go to 9000
664
end if
665
c
666
if (msglvl .gt. 2) then
667
call dvout (logfil, ncv, workl(iheigr), ndigit,
668
& '_neupd: Real part of the eigenvalues of H--reordered')
669
call dvout (logfil, ncv, workl(iheigi), ndigit,
670
& '_neupd: Imag part of the eigenvalues of H--reordered')
671
if (msglvl .gt. 3) then
672
call dmout (logfil, ncv, ncv, workl(iuptri), ldq,
673
& ndigit,
674
& '_neupd: Quasi-triangular matrix after re-ordering')
675
end if
676
end if
677
c
678
end if
679
c
680
c %---------------------------------------%
681
c | Copy the last row of the Schur vector |
682
c | into workl(ihbds). This will be used |
683
c | to compute the Ritz estimates of |
684
c | converged Ritz values. |
685
c %---------------------------------------%
686
c
687
call dcopy(ncv, workl(invsub+ncv-1), ldq, workl(ihbds), 1)
688
c
689
c %----------------------------------------------------%
690
c | Place the computed eigenvalues of H into DR and DI |
691
c | if a spectral transformation was not used. |
692
c %----------------------------------------------------%
693
c
694
if (type .eq. 'REGULR') then
695
call dcopy (nconv, workl(iheigr), 1, dr, 1)
696
call dcopy (nconv, workl(iheigi), 1, di, 1)
697
end if
698
c
699
c %----------------------------------------------------------%
700
c | Compute the QR factorization of the matrix representing |
701
c | the wanted invariant subspace located in the first NCONV |
702
c | columns of workl(invsub,ldq). |
703
c %----------------------------------------------------------%
704
c
705
call dgeqr2 (ncv, nconv, workl(invsub), ldq, workev,
706
& workev(ncv+1), ierr)
707
c
708
c %---------------------------------------------------------%
709
c | * Postmultiply V by Q using dorm2r. |
710
c | * Copy the first NCONV columns of VQ into Z. |
711
c | * Postmultiply Z by R. |
712
c | The N by NCONV matrix Z is now a matrix representation |
713
c | of the approximate invariant subspace associated with |
714
c | the Ritz values in workl(iheigr) and workl(iheigi) |
715
c | The first NCONV columns of V are now approximate Schur |
716
c | vectors associated with the real upper quasi-triangular |
717
c | matrix of order NCONV in workl(iuptri) |
718
c %---------------------------------------------------------%
719
c
720
call dorm2r ('Right', 'Notranspose', n, ncv, nconv,
721
& workl(invsub), ldq, workev, v, ldv, workd(n+1), ierr)
722
call dlacpy ('All', n, nconv, v, ldv, z, ldz)
723
c
724
do 20 j=1, nconv
725
c
726
c %---------------------------------------------------%
727
c | Perform both a column and row scaling if the |
728
c | diagonal element of workl(invsub,ldq) is negative |
729
c | I'm lazy and don't take advantage of the upper |
730
c | quasi-triangular form of workl(iuptri,ldq) |
731
c | Note that since Q is orthogonal, R is a diagonal |
732
c | matrix consisting of plus or minus ones |
733
c %---------------------------------------------------%
734
c
735
if (workl(invsub+(j-1)*ldq+j-1) .lt. zero) then
736
call dscal (nconv, -one, workl(iuptri+j-1), ldq)
737
call dscal (nconv, -one, workl(iuptri+(j-1)*ldq), 1)
738
end if
739
c
740
20 continue
741
c
742
if (howmny .eq. 'A') then
743
c
744
c %--------------------------------------------%
745
c | Compute the NCONV wanted eigenvectors of T |
746
c | located in workl(iuptri,ldq). |
747
c %--------------------------------------------%
748
c
749
do 30 j=1, ncv
750
if (j .le. nconv) then
751
select(j) = .true.
752
else
753
select(j) = .false.
754
end if
755
30 continue
756
c
757
call dtrevc ('Right', 'Select', select, ncv, workl(iuptri),
758
& ldq, vl, 1, workl(invsub), ldq, ncv, outncv, workev,
759
& ierr)
760
c
761
if (ierr .ne. 0) then
762
info = -9
763
go to 9000
764
end if
765
c
766
c %------------------------------------------------%
767
c | Scale the returning eigenvectors so that their |
768
c | Euclidean norms are all one. LAPACK subroutine |
769
c | dtrevc returns each eigenvector normalized so |
770
c | that the element of largest magnitude has |
771
c | magnitude 1; |
772
c %------------------------------------------------%
773
c
774
iconj = 0
775
do 40 j=1, nconv
776
c
777
if ( workl(iheigi+j-1) .eq. zero ) then
778
c
779
c %----------------------%
780
c | real eigenvalue case |
781
c %----------------------%
782
c
783
temp = dnrm2( ncv, workl(invsub+(j-1)*ldq), 1 )
784
call dscal ( ncv, one / temp,
785
& workl(invsub+(j-1)*ldq), 1 )
786
c
787
else
788
c
789
c %-------------------------------------------%
790
c | Complex conjugate pair case. Note that |
791
c | since the real and imaginary part of |
792
c | the eigenvector are stored in consecutive |
793
c | columns, we further normalize by the |
794
c | square root of two. |
795
c %-------------------------------------------%
796
c
797
if (iconj .eq. 0) then
798
temp = dlapy2( dnrm2( ncv, workl(invsub+(j-1)*ldq),
799
& 1 ), dnrm2( ncv, workl(invsub+j*ldq), 1) )
800
call dscal ( ncv, one / temp,
801
& workl(invsub+(j-1)*ldq), 1 )
802
call dscal ( ncv, one / temp,
803
& workl(invsub+j*ldq), 1 )
804
iconj = 1
805
else
806
iconj = 0
807
end if
808
c
809
end if
810
c
811
40 continue
812
c
813
call dgemv('T', ncv, nconv, one, workl(invsub),
814
& ldq, workl(ihbds), 1, zero, workev, 1)
815
c
816
iconj = 0
817
do 45 j=1, nconv
818
if (workl(iheigi+j-1) .ne. zero) then
819
c
820
c %-------------------------------------------%
821
c | Complex conjugate pair case. Note that |
822
c | since the real and imaginary part of |
823
c | the eigenvector are stored in consecutive |
824
c %-------------------------------------------%
825
c
826
if (iconj .eq. 0) then
827
workev(j) = dlapy2(workev(j), workev(j+1))
828
workev(j+1) = workev(j)
829
iconj = 1
830
else
831
iconj = 0
832
end if
833
end if
834
45 continue
835
c
836
if (msglvl .gt. 2) then
837
call dcopy(ncv, workl(invsub+ncv-1), ldq,
838
& workl(ihbds), 1)
839
call dvout (logfil, ncv, workl(ihbds), ndigit,
840
& '_neupd: Last row of the eigenvector matrix for T')
841
if (msglvl .gt. 3) then
842
call dmout (logfil, ncv, ncv, workl(invsub), ldq,
843
& ndigit, '_neupd: The eigenvector matrix for T')
844
end if
845
end if
846
c
847
c %---------------------------------------%
848
c | Copy Ritz estimates into workl(ihbds) |
849
c %---------------------------------------%
850
c
851
call dcopy(nconv, workev, 1, workl(ihbds), 1)
852
c
853
c %---------------------------------------------------------%
854
c | Compute the QR factorization of the eigenvector matrix |
855
c | associated with leading portion of T in the first NCONV |
856
c | columns of workl(invsub,ldq). |
857
c %---------------------------------------------------------%
858
c
859
call dgeqr2 (ncv, nconv, workl(invsub), ldq, workev,
860
& workev(ncv+1), ierr)
861
c
862
c %----------------------------------------------%
863
c | * Postmultiply Z by Q. |
864
c | * Postmultiply Z by R. |
865
c | The N by NCONV matrix Z is now contains the |
866
c | Ritz vectors associated with the Ritz values |
867
c | in workl(iheigr) and workl(iheigi). |
868
c %----------------------------------------------%
869
c
870
call dorm2r ('Right', 'Notranspose', n, ncv, nconv,
871
& workl(invsub), ldq, workev, z, ldz, workd(n+1), ierr)
872
c
873
call dtrmm ('Right', 'Upper', 'No transpose', 'Non-unit',
874
& n, nconv, one, workl(invsub), ldq, z, ldz)
875
c
876
end if
877
c
878
else
879
c
880
c %------------------------------------------------------%
881
c | An approximate invariant subspace is not needed. |
882
c | Place the Ritz values computed DNAUPD into DR and DI |
883
c %------------------------------------------------------%
884
c
885
call dcopy (nconv, workl(ritzr), 1, dr, 1)
886
call dcopy (nconv, workl(ritzi), 1, di, 1)
887
call dcopy (nconv, workl(ritzr), 1, workl(iheigr), 1)
888
call dcopy (nconv, workl(ritzi), 1, workl(iheigi), 1)
889
call dcopy (nconv, workl(bounds), 1, workl(ihbds), 1)
890
end if
891
c
892
c %------------------------------------------------%
893
c | Transform the Ritz values and possibly vectors |
894
c | and corresponding error bounds of OP to those |
895
c | of A*x = lambda*B*x. |
896
c %------------------------------------------------%
897
c
898
if (type .eq. 'REGULR') then
899
c
900
if (rvec)
901
& call dscal (ncv, rnorm, workl(ihbds), 1)
902
c
903
else
904
c
905
c %---------------------------------------%
906
c | A spectral transformation was used. |
907
c | * Determine the Ritz estimates of the |
908
c | Ritz values in the original system. |
909
c %---------------------------------------%
910
c
911
if (type .eq. 'SHIFTI') then
912
c
913
if (rvec)
914
& call dscal (ncv, rnorm, workl(ihbds), 1)
915
c
916
do 50 k=1, ncv
917
temp = dlapy2( workl(iheigr+k-1),
918
& workl(iheigi+k-1) )
919
workl(ihbds+k-1) = abs( workl(ihbds+k-1) )
920
& / temp / temp
921
50 continue
922
c
923
else if (type .eq. 'REALPT') then
924
c
925
do 60 k=1, ncv
926
60 continue
927
c
928
else if (type .eq. 'IMAGPT') then
929
c
930
do 70 k=1, ncv
931
70 continue
932
c
933
end if
934
c
935
c %-----------------------------------------------------------%
936
c | * Transform the Ritz values back to the original system. |
937
c | For TYPE = 'SHIFTI' the transformation is |
938
c | lambda = 1/theta + sigma |
939
c | For TYPE = 'REALPT' or 'IMAGPT' the user must from |
940
c | Rayleigh quotients or a projection. See remark 3 above.|
941
c | NOTES: |
942
c | *The Ritz vectors are not affected by the transformation. |
943
c %-----------------------------------------------------------%
944
c
945
if (type .eq. 'SHIFTI') then
946
c
947
do 80 k=1, ncv
948
temp = dlapy2( workl(iheigr+k-1),
949
& workl(iheigi+k-1) )
950
workl(iheigr+k-1) = workl(iheigr+k-1) / temp / temp
951
& + sigmar
952
workl(iheigi+k-1) = -workl(iheigi+k-1) / temp / temp
953
& + sigmai
954
80 continue
955
c
956
call dcopy (nconv, workl(iheigr), 1, dr, 1)
957
call dcopy (nconv, workl(iheigi), 1, di, 1)
958
c
959
else if (type .eq. 'REALPT' .or. type .eq. 'IMAGPT') then
960
c
961
call dcopy (nconv, workl(iheigr), 1, dr, 1)
962
call dcopy (nconv, workl(iheigi), 1, di, 1)
963
c
964
end if
965
c
966
end if
967
c
968
if (type .eq. 'SHIFTI' .and. msglvl .gt. 1) then
969
call dvout (logfil, nconv, dr, ndigit,
970
& '_neupd: Untransformed real part of the Ritz valuess.')
971
call dvout (logfil, nconv, di, ndigit,
972
& '_neupd: Untransformed imag part of the Ritz valuess.')
973
call dvout (logfil, nconv, workl(ihbds), ndigit,
974
& '_neupd: Ritz estimates of untransformed Ritz values.')
975
else if (type .eq. 'REGULR' .and. msglvl .gt. 1) then
976
call dvout (logfil, nconv, dr, ndigit,
977
& '_neupd: Real parts of converged Ritz values.')
978
call dvout (logfil, nconv, di, ndigit,
979
& '_neupd: Imag parts of converged Ritz values.')
980
call dvout (logfil, nconv, workl(ihbds), ndigit,
981
& '_neupd: Associated Ritz estimates.')
982
end if
983
c
984
c %-------------------------------------------------%
985
c | Eigenvector Purification step. Formally perform |
986
c | one of inverse subspace iteration. Only used |
987
c | for MODE = 2. |
988
c %-------------------------------------------------%
989
c
990
if (rvec .and. howmny .eq. 'A' .and. type .eq. 'SHIFTI') then
991
c
992
c %------------------------------------------------%
993
c | Purify the computed Ritz vectors by adding a |
994
c | little bit of the residual vector: |
995
c | T |
996
c | resid(:)*( e s ) / theta |
997
c | NCV |
998
c | where H s = s theta. Remember that when theta |
999
c | has nonzero imaginary part, the corresponding |
1000
c | Ritz vector is stored across two columns of Z. |
1001
c %------------------------------------------------%
1002
c
1003
iconj = 0
1004
do 110 j=1, nconv
1005
if (workl(iheigi+j-1) .eq. zero) then
1006
workev(j) = workl(invsub+(j-1)*ldq+ncv-1) /
1007
& workl(iheigr+j-1)
1008
else if (iconj .eq. 0) then
1009
temp = dlapy2( workl(iheigr+j-1), workl(iheigi+j-1) )
1010
workev(j) = ( workl(invsub+(j-1)*ldq+ncv-1) *
1011
& workl(iheigr+j-1) +
1012
& workl(invsub+j*ldq+ncv-1) *
1013
& workl(iheigi+j-1) ) / temp / temp
1014
workev(j+1) = ( workl(invsub+j*ldq+ncv-1) *
1015
& workl(iheigr+j-1) -
1016
& workl(invsub+(j-1)*ldq+ncv-1) *
1017
& workl(iheigi+j-1) ) / temp / temp
1018
iconj = 1
1019
else
1020
iconj = 0
1021
end if
1022
110 continue
1023
c
1024
c %---------------------------------------%
1025
c | Perform a rank one update to Z and |
1026
c | purify all the Ritz vectors together. |
1027
c %---------------------------------------%
1028
c
1029
call dger (n, nconv, one, resid, 1, workev, 1, z, ldz)
1030
c
1031
end if
1032
c
1033
9000 continue
1034
c
1035
return
1036
c
1037
c %---------------%
1038
c | End of DNEUPD |
1039
c %---------------%
1040
c
1041
end
1042
1043