Télécharger dneupd.eso

Retour à la liste

Numérotation des lignes :

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

© Cast3M 2003 - Tous droits réservés.
Mentions légales