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