You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

dgedmd.f90 46 kB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054
  1. SUBROUTINE DGEDMD( JOBS, JOBZ, JOBR, JOBF, WHTSVD, &
  2. M, N, X, LDX, Y, LDY, NRNK, TOL, &
  3. K, REIG, IMEIG, Z, LDZ, RES, &
  4. B, LDB, W, LDW, S, LDS, &
  5. WORK, LWORK, IWORK, LIWORK, INFO )
  6. ! March 2023
  7. !.....
  8. USE iso_fortran_env
  9. IMPLICIT NONE
  10. INTEGER, PARAMETER :: WP = real64
  11. !.....
  12. ! Scalar arguments
  13. CHARACTER, INTENT(IN) :: JOBS, JOBZ, JOBR, JOBF
  14. INTEGER, INTENT(IN) :: WHTSVD, M, N, LDX, LDY, &
  15. NRNK, LDZ, LDB, LDW, LDS, &
  16. LWORK, LIWORK
  17. INTEGER, INTENT(OUT) :: K, INFO
  18. REAL(KIND=WP), INTENT(IN) :: TOL
  19. ! Array arguments
  20. REAL(KIND=WP), INTENT(INOUT) :: X(LDX,*), Y(LDY,*)
  21. REAL(KIND=WP), INTENT(OUT) :: Z(LDZ,*), B(LDB,*), &
  22. W(LDW,*), S(LDS,*)
  23. REAL(KIND=WP), INTENT(OUT) :: REIG(*), IMEIG(*), &
  24. RES(*)
  25. REAL(KIND=WP), INTENT(OUT) :: WORK(*)
  26. INTEGER, INTENT(OUT) :: IWORK(*)
  27. !............................................................
  28. ! Purpose
  29. ! =======
  30. ! DGEDMD computes the Dynamic Mode Decomposition (DMD) for
  31. ! a pair of data snapshot matrices. For the input matrices
  32. ! X and Y such that Y = A*X with an unaccessible matrix
  33. ! A, DGEDMD computes a certain number of Ritz pairs of A using
  34. ! the standard Rayleigh-Ritz extraction from a subspace of
  35. ! range(X) that is determined using the leading left singular
  36. ! vectors of X. Optionally, DGEDMD returns the residuals
  37. ! of the computed Ritz pairs, the information needed for
  38. ! a refinement of the Ritz vectors, or the eigenvectors of
  39. ! the Exact DMD.
  40. ! For further details see the references listed
  41. ! below. For more details of the implementation see [3].
  42. !
  43. ! References
  44. ! ==========
  45. ! [1] P. Schmid: Dynamic mode decomposition of numerical
  46. ! and experimental data,
  47. ! Journal of Fluid Mechanics 656, 5-28, 2010.
  48. ! [2] Z. Drmac, I. Mezic, R. Mohr: Data driven modal
  49. ! decompositions: analysis and enhancements,
  50. ! SIAM J. on Sci. Comp. 40 (4), A2253-A2285, 2018.
  51. ! [3] Z. Drmac: A LAPACK implementation of the Dynamic
  52. ! Mode Decomposition I. Technical report. AIMDyn Inc.
  53. ! and LAPACK Working Note 298.
  54. ! [4] J. Tu, C. W. Rowley, D. M. Luchtenburg, S. L.
  55. ! Brunton, N. Kutz: On Dynamic Mode Decomposition:
  56. ! Theory and Applications, Journal of Computational
  57. ! Dynamics 1(2), 391 -421, 2014.
  58. !
  59. !......................................................................
  60. ! Developed and supported by:
  61. ! ===========================
  62. ! Developed and coded by Zlatko Drmac, Faculty of Science,
  63. ! University of Zagreb; drmac@math.hr
  64. ! In cooperation with
  65. ! AIMdyn Inc., Santa Barbara, CA.
  66. ! and supported by
  67. ! - DARPA SBIR project "Koopman Operator-Based Forecasting
  68. ! for Nonstationary Processes from Near-Term, Limited
  69. ! Observational Data" Contract No: W31P4Q-21-C-0007
  70. ! - DARPA PAI project "Physics-Informed Machine Learning
  71. ! Methodologies" Contract No: HR0011-18-9-0033
  72. ! - DARPA MoDyL project "A Data-Driven, Operator-Theoretic
  73. ! Framework for Space-Time Analysis of Process Dynamics"
  74. ! Contract No: HR0011-16-C-0116
  75. ! Any opinions, findings and conclusions or recommendations
  76. ! expressed in this material are those of the author and
  77. ! do not necessarily reflect the views of the DARPA SBIR
  78. ! Program Office
  79. !============================================================
  80. ! Distribution Statement A:
  81. ! Approved for Public Release, Distribution Unlimited.
  82. ! Cleared by DARPA on September 29, 2022
  83. !============================================================
  84. !............................................................
  85. ! Arguments
  86. ! =========
  87. ! JOBS (input) CHARACTER*1
  88. ! Determines whether the initial data snapshots are scaled
  89. ! by a diagonal matrix.
  90. ! 'S' :: The data snapshots matrices X and Y are multiplied
  91. ! with a diagonal matrix D so that X*D has unit
  92. ! nonzero columns (in the Euclidean 2-norm)
  93. ! 'C' :: The snapshots are scaled as with the 'S' option.
  94. ! If it is found that an i-th column of X is zero
  95. ! vector and the corresponding i-th column of Y is
  96. ! non-zero, then the i-th column of Y is set to
  97. ! zero and a warning flag is raised.
  98. ! 'Y' :: The data snapshots matrices X and Y are multiplied
  99. ! by a diagonal matrix D so that Y*D has unit
  100. ! nonzero columns (in the Euclidean 2-norm)
  101. ! 'N' :: No data scaling.
  102. !.....
  103. ! JOBZ (input) CHARACTER*1
  104. ! Determines whether the eigenvectors (Koopman modes) will
  105. ! be computed.
  106. ! 'V' :: The eigenvectors (Koopman modes) will be computed
  107. ! and returned in the matrix Z.
  108. ! See the description of Z.
  109. ! 'F' :: The eigenvectors (Koopman modes) will be returned
  110. ! in factored form as the product X(:,1:K)*W, where X
  111. ! contains a POD basis (leading left singular vectors
  112. ! of the data matrix X) and W contains the eigenvectors
  113. ! of the corresponding Rayleigh quotient.
  114. ! See the descriptions of K, X, W, Z.
  115. ! 'N' :: The eigenvectors are not computed.
  116. !.....
  117. ! JOBR (input) CHARACTER*1
  118. ! Determines whether to compute the residuals.
  119. ! 'R' :: The residuals for the computed eigenpairs will be
  120. ! computed and stored in the array RES.
  121. ! See the description of RES.
  122. ! For this option to be legal, JOBZ must be 'V'.
  123. ! 'N' :: The residuals are not computed.
  124. !.....
  125. ! JOBF (input) CHARACTER*1
  126. ! Specifies whether to store information needed for post-
  127. ! processing (e.g. computing refined Ritz vectors)
  128. ! 'R' :: The matrix needed for the refinement of the Ritz
  129. ! vectors is computed and stored in the array B.
  130. ! See the description of B.
  131. ! 'E' :: The unscaled eigenvectors of the Exact DMD are
  132. ! computed and returned in the array B. See the
  133. ! description of B.
  134. ! 'N' :: No eigenvector refinement data is computed.
  135. !.....
  136. ! WHTSVD (input) INTEGER, WHSTVD in { 1, 2, 3, 4 }
  137. ! Allows for a selection of the SVD algorithm from the
  138. ! LAPACK library.
  139. ! 1 :: DGESVD (the QR SVD algorithm)
  140. ! 2 :: DGESDD (the Divide and Conquer algorithm; if enough
  141. ! workspace available, this is the fastest option)
  142. ! 3 :: DGESVDQ (the preconditioned QR SVD ; this and 4
  143. ! are the most accurate options)
  144. ! 4 :: DGEJSV (the preconditioned Jacobi SVD; this and 3
  145. ! are the most accurate options)
  146. ! For the four methods above, a significant difference in
  147. ! the accuracy of small singular values is possible if
  148. ! the snapshots vary in norm so that X is severely
  149. ! ill-conditioned. If small (smaller than EPS*||X||)
  150. ! singular values are of interest and JOBS=='N', then
  151. ! the options (3, 4) give the most accurate results, where
  152. ! the option 4 is slightly better and with stronger
  153. ! theoretical background.
  154. ! If JOBS=='S', i.e. the columns of X will be normalized,
  155. ! then all methods give nearly equally accurate results.
  156. !.....
  157. ! M (input) INTEGER, M>= 0
  158. ! The state space dimension (the row dimension of X, Y).
  159. !.....
  160. ! N (input) INTEGER, 0 <= N <= M
  161. ! The number of data snapshot pairs
  162. ! (the number of columns of X and Y).
  163. !.....
  164. ! X (input/output) REAL(KIND=WP) M-by-N array
  165. ! > On entry, X contains the data snapshot matrix X. It is
  166. ! assumed that the column norms of X are in the range of
  167. ! the normalized floating point numbers.
  168. ! < On exit, the leading K columns of X contain a POD basis,
  169. ! i.e. the leading K left singular vectors of the input
  170. ! data matrix X, U(:,1:K). All N columns of X contain all
  171. ! left singular vectors of the input matrix X.
  172. ! See the descriptions of K, Z and W.
  173. !.....
  174. ! LDX (input) INTEGER, LDX >= M
  175. ! The leading dimension of the array X.
  176. !.....
  177. ! Y (input/workspace/output) REAL(KIND=WP) M-by-N array
  178. ! > On entry, Y contains the data snapshot matrix Y
  179. ! < On exit,
  180. ! If JOBR == 'R', the leading K columns of Y contain
  181. ! the residual vectors for the computed Ritz pairs.
  182. ! See the description of RES.
  183. ! If JOBR == 'N', Y contains the original input data,
  184. ! scaled according to the value of JOBS.
  185. !.....
  186. ! LDY (input) INTEGER , LDY >= M
  187. ! The leading dimension of the array Y.
  188. !.....
  189. ! NRNK (input) INTEGER
  190. ! Determines the mode how to compute the numerical rank,
  191. ! i.e. how to truncate small singular values of the input
  192. ! matrix X. On input, if
  193. ! NRNK = -1 :: i-th singular value sigma(i) is truncated
  194. ! if sigma(i) <= TOL*sigma(1).
  195. ! This option is recommended.
  196. ! NRNK = -2 :: i-th singular value sigma(i) is truncated
  197. ! if sigma(i) <= TOL*sigma(i-1)
  198. ! This option is included for R&D purposes.
  199. ! It requires highly accurate SVD, which
  200. ! may not be feasible.
  201. !
  202. ! The numerical rank can be enforced by using positive
  203. ! value of NRNK as follows:
  204. ! 0 < NRNK <= N :: at most NRNK largest singular values
  205. ! will be used. If the number of the computed nonzero
  206. ! singular values is less than NRNK, then only those
  207. ! nonzero values will be used and the actually used
  208. ! dimension is less than NRNK. The actual number of
  209. ! the nonzero singular values is returned in the variable
  210. ! K. See the descriptions of TOL and K.
  211. !.....
  212. ! TOL (input) REAL(KIND=WP), 0 <= TOL < 1
  213. ! The tolerance for truncating small singular values.
  214. ! See the description of NRNK.
  215. !.....
  216. ! K (output) INTEGER, 0 <= K <= N
  217. ! The dimension of the POD basis for the data snapshot
  218. ! matrix X and the number of the computed Ritz pairs.
  219. ! The value of K is determined according to the rule set
  220. ! by the parameters NRNK and TOL.
  221. ! See the descriptions of NRNK and TOL.
  222. !.....
  223. ! REIG (output) REAL(KIND=WP) N-by-1 array
  224. ! The leading K (K<=N) entries of REIG contain
  225. ! the real parts of the computed eigenvalues
  226. ! REIG(1:K) + sqrt(-1)*IMEIG(1:K).
  227. ! See the descriptions of K, IMEIG, and Z.
  228. !.....
  229. ! IMEIG (output) REAL(KIND=WP) N-by-1 array
  230. ! The leading K (K<=N) entries of IMEIG contain
  231. ! the imaginary parts of the computed eigenvalues
  232. ! REIG(1:K) + sqrt(-1)*IMEIG(1:K).
  233. ! The eigenvalues are determined as follows:
  234. ! If IMEIG(i) == 0, then the corresponding eigenvalue is
  235. ! real, LAMBDA(i) = REIG(i).
  236. ! If IMEIG(i)>0, then the corresponding complex
  237. ! conjugate pair of eigenvalues reads
  238. ! LAMBDA(i) = REIG(i) + sqrt(-1)*IMAG(i)
  239. ! LAMBDA(i+1) = REIG(i) - sqrt(-1)*IMAG(i)
  240. ! That is, complex conjugate pairs have consecutive
  241. ! indices (i,i+1), with the positive imaginary part
  242. ! listed first.
  243. ! See the descriptions of K, REIG, and Z.
  244. !.....
  245. ! Z (workspace/output) REAL(KIND=WP) M-by-N array
  246. ! If JOBZ =='V' then
  247. ! Z contains real Ritz vectors as follows:
  248. ! If IMEIG(i)=0, then Z(:,i) is an eigenvector of
  249. ! the i-th Ritz value; ||Z(:,i)||_2=1.
  250. ! If IMEIG(i) > 0 (and IMEIG(i+1) < 0) then
  251. ! [Z(:,i) Z(:,i+1)] span an invariant subspace and
  252. ! the Ritz values extracted from this subspace are
  253. ! REIG(i) + sqrt(-1)*IMEIG(i) and
  254. ! REIG(i) - sqrt(-1)*IMEIG(i).
  255. ! The corresponding eigenvectors are
  256. ! Z(:,i) + sqrt(-1)*Z(:,i+1) and
  257. ! Z(:,i) - sqrt(-1)*Z(:,i+1), respectively.
  258. ! || Z(:,i:i+1)||_F = 1.
  259. ! If JOBZ == 'F', then the above descriptions hold for
  260. ! the columns of X(:,1:K)*W(1:K,1:K), where the columns
  261. ! of W(1:k,1:K) are the computed eigenvectors of the
  262. ! K-by-K Rayleigh quotient. The columns of W(1:K,1:K)
  263. ! are similarly structured: If IMEIG(i) == 0 then
  264. ! X(:,1:K)*W(:,i) is an eigenvector, and if IMEIG(i)>0
  265. ! then X(:,1:K)*W(:,i)+sqrt(-1)*X(:,1:K)*W(:,i+1) and
  266. ! X(:,1:K)*W(:,i)-sqrt(-1)*X(:,1:K)*W(:,i+1)
  267. ! are the eigenvectors of LAMBDA(i), LAMBDA(i+1).
  268. ! See the descriptions of REIG, IMEIG, X and W.
  269. !.....
  270. ! LDZ (input) INTEGER , LDZ >= M
  271. ! The leading dimension of the array Z.
  272. !.....
  273. ! RES (output) REAL(KIND=WP) N-by-1 array
  274. ! RES(1:K) contains the residuals for the K computed
  275. ! Ritz pairs.
  276. ! If LAMBDA(i) is real, then
  277. ! RES(i) = || A * Z(:,i) - LAMBDA(i)*Z(:,i))||_2.
  278. ! If [LAMBDA(i), LAMBDA(i+1)] is a complex conjugate pair
  279. ! then
  280. ! RES(i)=RES(i+1) = || A * Z(:,i:i+1) - Z(:,i:i+1) *B||_F
  281. ! where B = [ real(LAMBDA(i)) imag(LAMBDA(i)) ]
  282. ! [-imag(LAMBDA(i)) real(LAMBDA(i)) ].
  283. ! It holds that
  284. ! RES(i) = || A*ZC(:,i) - LAMBDA(i) *ZC(:,i) ||_2
  285. ! RES(i+1) = || A*ZC(:,i+1) - LAMBDA(i+1)*ZC(:,i+1) ||_2
  286. ! where ZC(:,i) = Z(:,i) + sqrt(-1)*Z(:,i+1)
  287. ! ZC(:,i+1) = Z(:,i) - sqrt(-1)*Z(:,i+1)
  288. ! See the description of REIG, IMEIG and Z.
  289. !.....
  290. ! B (output) REAL(KIND=WP) M-by-N array.
  291. ! IF JOBF =='R', B(1:M,1:K) contains A*U(:,1:K), and can
  292. ! be used for computing the refined vectors; see further
  293. ! details in the provided references.
  294. ! If JOBF == 'E', B(1:M,1;K) contains
  295. ! A*U(:,1:K)*W(1:K,1:K), which are the vectors from the
  296. ! Exact DMD, up to scaling by the inverse eigenvalues.
  297. ! If JOBF =='N', then B is not referenced.
  298. ! See the descriptions of X, W, K.
  299. !.....
  300. ! LDB (input) INTEGER, LDB >= M
  301. ! The leading dimension of the array B.
  302. !.....
  303. ! W (workspace/output) REAL(KIND=WP) N-by-N array
  304. ! On exit, W(1:K,1:K) contains the K computed
  305. ! eigenvectors of the matrix Rayleigh quotient (real and
  306. ! imaginary parts for each complex conjugate pair of the
  307. ! eigenvalues). The Ritz vectors (returned in Z) are the
  308. ! product of X (containing a POD basis for the input
  309. ! matrix X) and W. See the descriptions of K, S, X and Z.
  310. ! W is also used as a workspace to temporarily store the
  311. ! right singular vectors of X.
  312. !.....
  313. ! LDW (input) INTEGER, LDW >= N
  314. ! The leading dimension of the array W.
  315. !.....
  316. ! S (workspace/output) REAL(KIND=WP) N-by-N array
  317. ! The array S(1:K,1:K) is used for the matrix Rayleigh
  318. ! quotient. This content is overwritten during
  319. ! the eigenvalue decomposition by DGEEV.
  320. ! See the description of K.
  321. !.....
  322. ! LDS (input) INTEGER, LDS >= N
  323. ! The leading dimension of the array S.
  324. !.....
  325. ! WORK (workspace/output) REAL(KIND=WP) LWORK-by-1 array
  326. ! On exit, WORK(1:N) contains the singular values of
  327. ! X (for JOBS=='N') or column scaled X (JOBS=='S', 'C').
  328. ! If WHTSVD==4, then WORK(N+1) and WORK(N+2) contain
  329. ! scaling factor WORK(N+2)/WORK(N+1) used to scale X
  330. ! and Y to avoid overflow in the SVD of X.
  331. ! This may be of interest if the scaling option is off
  332. ! and as many as possible smallest eigenvalues are
  333. ! desired to the highest feasible accuracy.
  334. ! If the call to DGEDMD is only workspace query, then
  335. ! WORK(1) contains the minimal workspace length and
  336. ! WORK(2) is the optimal workspace length. Hence, the
  337. ! leng of work is at least 2.
  338. ! See the description of LWORK.
  339. !.....
  340. ! LWORK (input) INTEGER
  341. ! The minimal length of the workspace vector WORK.
  342. ! LWORK is calculated as follows:
  343. ! If WHTSVD == 1 ::
  344. ! If JOBZ == 'V', then
  345. ! LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,4*N)).
  346. ! If JOBZ == 'N' then
  347. ! LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,3*N)).
  348. ! Here LWORK_SVD = MAX(1,3*N+M,5*N) is the minimal
  349. ! workspace length of DGESVD.
  350. ! If WHTSVD == 2 ::
  351. ! If JOBZ == 'V', then
  352. ! LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,4*N))
  353. ! If JOBZ == 'N', then
  354. ! LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,3*N))
  355. ! Here LWORK_SVD = MAX(M, 5*N*N+4*N)+3*N*N is the
  356. ! minimal workspace length of DGESDD.
  357. ! If WHTSVD == 3 ::
  358. ! If JOBZ == 'V', then
  359. ! LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,4*N))
  360. ! If JOBZ == 'N', then
  361. ! LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,3*N))
  362. ! Here LWORK_SVD = N+M+MAX(3*N+1,
  363. ! MAX(1,3*N+M,5*N),MAX(1,N))
  364. ! is the minimal workspace length of DGESVDQ.
  365. ! If WHTSVD == 4 ::
  366. ! If JOBZ == 'V', then
  367. ! LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,4*N))
  368. ! If JOBZ == 'N', then
  369. ! LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,3*N))
  370. ! Here LWORK_SVD = MAX(7,2*M+N,6*N+2*N*N) is the
  371. ! minimal workspace length of DGEJSV.
  372. ! The above expressions are not simplified in order to
  373. ! make the usage of WORK more transparent, and for
  374. ! easier checking. In any case, LWORK >= 2.
  375. ! If on entry LWORK = -1, then a workspace query is
  376. ! assumed and the procedure only computes the minimal
  377. ! and the optimal workspace lengths for both WORK and
  378. ! IWORK. See the descriptions of WORK and IWORK.
  379. !.....
  380. ! IWORK (workspace/output) INTEGER LIWORK-by-1 array
  381. ! Workspace that is required only if WHTSVD equals
  382. ! 2 , 3 or 4. (See the description of WHTSVD).
  383. ! If on entry LWORK =-1 or LIWORK=-1, then the
  384. ! minimal length of IWORK is computed and returned in
  385. ! IWORK(1). See the description of LIWORK.
  386. !.....
  387. ! LIWORK (input) INTEGER
  388. ! The minimal length of the workspace vector IWORK.
  389. ! If WHTSVD == 1, then only IWORK(1) is used; LIWORK >=1
  390. ! If WHTSVD == 2, then LIWORK >= MAX(1,8*MIN(M,N))
  391. ! If WHTSVD == 3, then LIWORK >= MAX(1,M+N-1)
  392. ! If WHTSVD == 4, then LIWORK >= MAX(3,M+3*N)
  393. ! If on entry LIWORK = -1, then a workspace query is
  394. ! assumed and the procedure only computes the minimal
  395. ! and the optimal workspace lengths for both WORK and
  396. ! IWORK. See the descriptions of WORK and IWORK.
  397. !.....
  398. ! INFO (output) INTEGER
  399. ! -i < 0 :: On entry, the i-th argument had an
  400. ! illegal value
  401. ! = 0 :: Successful return.
  402. ! = 1 :: Void input. Quick exit (M=0 or N=0).
  403. ! = 2 :: The SVD computation of X did not converge.
  404. ! Suggestion: Check the input data and/or
  405. ! repeat with different WHTSVD.
  406. ! = 3 :: The computation of the eigenvalues did not
  407. ! converge.
  408. ! = 4 :: If data scaling was requested on input and
  409. ! the procedure found inconsistency in the data
  410. ! such that for some column index i,
  411. ! X(:,i) = 0 but Y(:,i) /= 0, then Y(:,i) is set
  412. ! to zero if JOBS=='C'. The computation proceeds
  413. ! with original or modified data and warning
  414. ! flag is set with INFO=4.
  415. !.............................................................
  416. !.............................................................
  417. ! Parameters
  418. ! ~~~~~~~~~~
  419. REAL(KIND=WP), PARAMETER :: ONE = 1.0_WP
  420. REAL(KIND=WP), PARAMETER :: ZERO = 0.0_WP
  421. ! Local scalars
  422. ! ~~~~~~~~~~~~~
  423. REAL(KIND=WP) :: OFL, ROOTSC, SCALE, SMALL, &
  424. SSUM, XSCL1, XSCL2
  425. INTEGER :: i, j, IMINWR, INFO1, INFO2, &
  426. LWRKEV, LWRSDD, LWRSVD, &
  427. LWRSVQ, MLWORK, MWRKEV, MWRSDD, &
  428. MWRSVD, MWRSVJ, MWRSVQ, NUMRNK, &
  429. OLWORK
  430. LOGICAL :: BADXY, LQUERY, SCCOLX, SCCOLY, &
  431. WNTEX, WNTREF, WNTRES, WNTVEC
  432. CHARACTER :: JOBZL, T_OR_N
  433. CHARACTER :: JSVOPT
  434. ! Local arrays
  435. ! ~~~~~~~~~~~~
  436. REAL(KIND=WP) :: AB(2,2), RDUMMY(2), RDUMMY2(2)
  437. ! External functions (BLAS and LAPACK)
  438. ! ~~~~~~~~~~~~~~~~~
  439. REAL(KIND=WP) DLANGE, DLAMCH, DNRM2
  440. EXTERNAL DLANGE, DLAMCH, DNRM2, IDAMAX
  441. INTEGER IDAMAX
  442. LOGICAL DISNAN, LSAME
  443. EXTERNAL DISNAN, LSAME
  444. ! External subroutines (BLAS and LAPACK)
  445. ! ~~~~~~~~~~~~~~~~~~~~
  446. EXTERNAL DAXPY, DGEMM, DSCAL
  447. EXTERNAL DGEEV, DGEJSV, DGESDD, DGESVD, DGESVDQ, &
  448. DLACPY, DLASCL, DLASSQ, XERBLA
  449. ! Intrinsic functions
  450. ! ~~~~~~~~~~~~~~~~~~~
  451. INTRINSIC DBLE, INT, MAX, SQRT
  452. !............................................................
  453. !
  454. ! Test the input arguments
  455. !
  456. WNTRES = LSAME(JOBR,'R')
  457. SCCOLX = LSAME(JOBS,'S') .OR. LSAME(JOBS,'C')
  458. SCCOLY = LSAME(JOBS,'Y')
  459. WNTVEC = LSAME(JOBZ,'V')
  460. WNTREF = LSAME(JOBF,'R')
  461. WNTEX = LSAME(JOBF,'E')
  462. INFO = 0
  463. LQUERY = ( ( LWORK == -1 ) .OR. ( LIWORK == -1 ) )
  464. !
  465. IF ( .NOT. (SCCOLX .OR. SCCOLY .OR. &
  466. LSAME(JOBS,'N')) ) THEN
  467. INFO = -1
  468. ELSE IF ( .NOT. (WNTVEC .OR. LSAME(JOBZ,'N') &
  469. .OR. LSAME(JOBZ,'F')) ) THEN
  470. INFO = -2
  471. ELSE IF ( .NOT. (WNTRES .OR. LSAME(JOBR,'N')) .OR. &
  472. ( WNTRES .AND. (.NOT.WNTVEC) ) ) THEN
  473. INFO = -3
  474. ELSE IF ( .NOT. (WNTREF .OR. WNTEX .OR. &
  475. LSAME(JOBF,'N') ) ) THEN
  476. INFO = -4
  477. ELSE IF ( .NOT.((WHTSVD == 1) .OR. (WHTSVD == 2) .OR. &
  478. (WHTSVD == 3) .OR. (WHTSVD == 4) )) THEN
  479. INFO = -5
  480. ELSE IF ( M < 0 ) THEN
  481. INFO = -6
  482. ELSE IF ( ( N < 0 ) .OR. ( N > M ) ) THEN
  483. INFO = -7
  484. ELSE IF ( LDX < M ) THEN
  485. INFO = -9
  486. ELSE IF ( LDY < M ) THEN
  487. INFO = -11
  488. ELSE IF ( .NOT. (( NRNK == -2).OR.(NRNK == -1).OR. &
  489. ((NRNK >= 1).AND.(NRNK <=N ))) ) THEN
  490. INFO = -12
  491. ELSE IF ( ( TOL < ZERO ) .OR. ( TOL >= ONE ) ) THEN
  492. INFO = -13
  493. ELSE IF ( LDZ < M ) THEN
  494. INFO = -18
  495. ELSE IF ( (WNTREF .OR. WNTEX ) .AND. ( LDB < M ) ) THEN
  496. INFO = -21
  497. ELSE IF ( LDW < N ) THEN
  498. INFO = -23
  499. ELSE IF ( LDS < N ) THEN
  500. INFO = -25
  501. END IF
  502. !
  503. IF ( INFO == 0 ) THEN
  504. ! Compute the minimal and the optimal workspace
  505. ! requirements. Simulate running the code and
  506. ! determine minimal and optimal sizes of the
  507. ! workspace at any moment of the run.
  508. IF ( N == 0 ) THEN
  509. ! Quick return. All output except K is void.
  510. ! INFO=1 signals the void input.
  511. ! In case of a workspace query, the default
  512. ! minimal workspace lengths are returned.
  513. IF ( LQUERY ) THEN
  514. IWORK(1) = 1
  515. WORK(1) = 2
  516. WORK(2) = 2
  517. ELSE
  518. K = 0
  519. END IF
  520. INFO = 1
  521. RETURN
  522. END IF
  523. MLWORK = MAX(2,N)
  524. OLWORK = MAX(2,N)
  525. IMINWR = 1
  526. SELECT CASE ( WHTSVD )
  527. CASE (1)
  528. ! The following is specified as the minimal
  529. ! length of WORK in the definition of DGESVD:
  530. ! MWRSVD = MAX(1,3*MIN(M,N)+MAX(M,N),5*MIN(M,N))
  531. MWRSVD = MAX(1,3*MIN(M,N)+MAX(M,N),5*MIN(M,N))
  532. MLWORK = MAX(MLWORK,N + MWRSVD)
  533. IF ( LQUERY ) THEN
  534. CALL DGESVD( 'O', 'S', M, N, X, LDX, WORK, &
  535. B, LDB, W, LDW, RDUMMY, -1, INFO1 )
  536. LWRSVD = MAX( MWRSVD, INT( RDUMMY(1) ) )
  537. OLWORK = MAX(OLWORK,N + LWRSVD)
  538. END IF
  539. CASE (2)
  540. ! The following is specified as the minimal
  541. ! length of WORK in the definition of DGESDD:
  542. ! MWRSDD = 3*MIN(M,N)*MIN(M,N) +
  543. ! MAX( MAX(M,N),5*MIN(M,N)*MIN(M,N)+4*MIN(M,N) )
  544. ! IMINWR = 8*MIN(M,N)
  545. MWRSDD = 3*MIN(M,N)*MIN(M,N) + &
  546. MAX( MAX(M,N),5*MIN(M,N)*MIN(M,N)+4*MIN(M,N) )
  547. MLWORK = MAX(MLWORK,N + MWRSDD)
  548. IMINWR = 8*MIN(M,N)
  549. IF ( LQUERY ) THEN
  550. CALL DGESDD( 'O', M, N, X, LDX, WORK, B, &
  551. LDB, W, LDW, RDUMMY, -1, IWORK, INFO1 )
  552. LWRSDD = MAX( MWRSDD, INT( RDUMMY(1) ) )
  553. OLWORK = MAX(OLWORK,N + LWRSDD)
  554. END IF
  555. CASE (3)
  556. !LWQP3 = 3*N+1
  557. !LWORQ = MAX(N, 1)
  558. !MWRSVD = MAX(1,3*MIN(M,N)+MAX(M,N),5*MIN(M,N))
  559. !MWRSVQ = N + MAX( LWQP3, MWRSVD, LWORQ ) + MAX(M,2)
  560. !MLWORK = N + MWRSVQ
  561. !IMINWR = M+N-1
  562. CALL DGESVDQ( 'H', 'P', 'N', 'R', 'R', M, N, &
  563. X, LDX, WORK, Z, LDZ, W, LDW, &
  564. NUMRNK, IWORK, LIWORK, RDUMMY, &
  565. -1, RDUMMY2, -1, INFO1 )
  566. IMINWR = IWORK(1)
  567. MWRSVQ = INT(RDUMMY(2))
  568. MLWORK = MAX(MLWORK,N+MWRSVQ+INT(RDUMMY2(1)))
  569. IF ( LQUERY ) THEN
  570. LWRSVQ = MAX( MWRSVQ, INT(RDUMMY(1)) )
  571. OLWORK = MAX(OLWORK,N+LWRSVQ+INT(RDUMMY2(1)))
  572. END IF
  573. CASE (4)
  574. JSVOPT = 'J'
  575. !MWRSVJ = MAX( 7, 2*M+N, 6*N+2*N*N ) ! for JSVOPT='V'
  576. MWRSVJ = MAX( 7, 2*M+N, 4*N+N*N, 2*N+N*N+6 )
  577. MLWORK = MAX(MLWORK,N+MWRSVJ)
  578. IMINWR = MAX( 3, M+3*N )
  579. IF ( LQUERY ) THEN
  580. OLWORK = MAX(OLWORK,N+MWRSVJ)
  581. END IF
  582. END SELECT
  583. IF ( WNTVEC .OR. WNTEX .OR. LSAME(JOBZ,'F') ) THEN
  584. JOBZL = 'V'
  585. ELSE
  586. JOBZL = 'N'
  587. END IF
  588. ! Workspace calculation to the DGEEV call
  589. IF ( LSAME(JOBZL,'V') ) THEN
  590. MWRKEV = MAX( 1, 4*N )
  591. ELSE
  592. MWRKEV = MAX( 1, 3*N )
  593. END IF
  594. MLWORK = MAX(MLWORK,N+MWRKEV)
  595. IF ( LQUERY ) THEN
  596. CALL DGEEV( 'N', JOBZL, N, S, LDS, REIG, &
  597. IMEIG, W, LDW, W, LDW, RDUMMY, -1, INFO1 )
  598. LWRKEV = MAX( MWRKEV, INT(RDUMMY(1)) )
  599. OLWORK = MAX( OLWORK, N+LWRKEV )
  600. END IF
  601. !
  602. IF ( LIWORK < IMINWR .AND. (.NOT.LQUERY) ) INFO = -29
  603. IF ( LWORK < MLWORK .AND. (.NOT.LQUERY) ) INFO = -27
  604. END IF
  605. !
  606. IF( INFO /= 0 ) THEN
  607. CALL XERBLA( 'DGEDMD', -INFO )
  608. RETURN
  609. ELSE IF ( LQUERY ) THEN
  610. ! Return minimal and optimal workspace sizes
  611. IWORK(1) = IMINWR
  612. WORK(1) = MLWORK
  613. WORK(2) = OLWORK
  614. RETURN
  615. END IF
  616. !............................................................
  617. !
  618. OFL = DLAMCH('O')
  619. SMALL = DLAMCH('S')
  620. BADXY = .FALSE.
  621. !
  622. ! <1> Optional scaling of the snapshots (columns of X, Y)
  623. ! ==========================================================
  624. IF ( SCCOLX ) THEN
  625. ! The columns of X will be normalized.
  626. ! To prevent overflows, the column norms of X are
  627. ! carefully computed using DLASSQ.
  628. K = 0
  629. DO i = 1, N
  630. !WORK(i) = DNRM2( M, X(1,i), 1 )
  631. SCALE = ZERO
  632. CALL DLASSQ( M, X(1,i), 1, SCALE, SSUM )
  633. IF ( DISNAN(SCALE) .OR. DISNAN(SSUM) ) THEN
  634. K = 0
  635. INFO = -8
  636. CALL XERBLA('DGEDMD',-INFO)
  637. END IF
  638. IF ( (SCALE /= ZERO) .AND. (SSUM /= ZERO) ) THEN
  639. ROOTSC = SQRT(SSUM)
  640. IF ( SCALE .GE. (OFL / ROOTSC) ) THEN
  641. ! Norm of X(:,i) overflows. First, X(:,i)
  642. ! is scaled by
  643. ! ( ONE / ROOTSC ) / SCALE = 1/||X(:,i)||_2.
  644. ! Next, the norm of X(:,i) is stored without
  645. ! overflow as WORK(i) = - SCALE * (ROOTSC/M),
  646. ! the minus sign indicating the 1/M factor.
  647. ! Scaling is performed without overflow, and
  648. ! underflow may occur in the smallest entries
  649. ! of X(:,i). The relative backward and forward
  650. ! errors are small in the ell_2 norm.
  651. CALL DLASCL( 'G', 0, 0, SCALE, ONE/ROOTSC, &
  652. M, 1, X(1,i), M, INFO2 )
  653. WORK(i) = - SCALE * ( ROOTSC / DBLE(M) )
  654. ELSE
  655. ! X(:,i) will be scaled to unit 2-norm
  656. WORK(i) = SCALE * ROOTSC
  657. CALL DLASCL( 'G',0, 0, WORK(i), ONE, M, 1, &
  658. X(1,i), M, INFO2 ) ! LAPACK CALL
  659. ! X(1:M,i) = (ONE/WORK(i)) * X(1:M,i) ! INTRINSIC
  660. END IF
  661. ELSE
  662. WORK(i) = ZERO
  663. K = K + 1
  664. END IF
  665. END DO
  666. IF ( K == N ) THEN
  667. ! All columns of X are zero. Return error code -8.
  668. ! (the 8th input variable had an illegal value)
  669. K = 0
  670. INFO = -8
  671. CALL XERBLA('DGEDMD',-INFO)
  672. RETURN
  673. END IF
  674. DO i = 1, N
  675. ! Now, apply the same scaling to the columns of Y.
  676. IF ( WORK(i) > ZERO ) THEN
  677. CALL DSCAL( M, ONE/WORK(i), Y(1,i), 1 ) ! BLAS CALL
  678. ! Y(1:M,i) = (ONE/WORK(i)) * Y(1:M,i) ! INTRINSIC
  679. ELSE IF ( WORK(i) < ZERO ) THEN
  680. CALL DLASCL( 'G', 0, 0, -WORK(i), &
  681. ONE/DBLE(M), M, 1, Y(1,i), M, INFO2 ) ! LAPACK CALL
  682. ELSE IF ( Y(IDAMAX(M, Y(1,i),1),i ) &
  683. /= ZERO ) THEN
  684. ! X(:,i) is zero vector. For consistency,
  685. ! Y(:,i) should also be zero. If Y(:,i) is not
  686. ! zero, then the data might be inconsistent or
  687. ! corrupted. If JOBS == 'C', Y(:,i) is set to
  688. ! zero and a warning flag is raised.
  689. ! The computation continues but the
  690. ! situation will be reported in the output.
  691. BADXY = .TRUE.
  692. IF ( LSAME(JOBS,'C')) &
  693. CALL DSCAL( M, ZERO, Y(1,i), 1 ) ! BLAS CALL
  694. END IF
  695. END DO
  696. END IF
  697. !
  698. IF ( SCCOLY ) THEN
  699. ! The columns of Y will be normalized.
  700. ! To prevent overflows, the column norms of Y are
  701. ! carefully computed using DLASSQ.
  702. DO i = 1, N
  703. !WORK(i) = DNRM2( M, Y(1,i), 1 )
  704. SCALE = ZERO
  705. CALL DLASSQ( M, Y(1,i), 1, SCALE, SSUM )
  706. IF ( DISNAN(SCALE) .OR. DISNAN(SSUM) ) THEN
  707. K = 0
  708. INFO = -10
  709. CALL XERBLA('DGEDMD',-INFO)
  710. END IF
  711. IF ( SCALE /= ZERO .AND. (SSUM /= ZERO) ) THEN
  712. ROOTSC = SQRT(SSUM)
  713. IF ( SCALE .GE. (OFL / ROOTSC) ) THEN
  714. ! Norm of Y(:,i) overflows. First, Y(:,i)
  715. ! is scaled by
  716. ! ( ONE / ROOTSC ) / SCALE = 1/||Y(:,i)||_2.
  717. ! Next, the norm of Y(:,i) is stored without
  718. ! overflow as WORK(i) = - SCALE * (ROOTSC/M),
  719. ! the minus sign indicating the 1/M factor.
  720. ! Scaling is performed without overflow, and
  721. ! underflow may occur in the smallest entries
  722. ! of Y(:,i). The relative backward and forward
  723. ! errors are small in the ell_2 norm.
  724. CALL DLASCL( 'G', 0, 0, SCALE, ONE/ROOTSC, &
  725. M, 1, Y(1,i), M, INFO2 )
  726. WORK(i) = - SCALE * ( ROOTSC / DBLE(M) )
  727. ELSE
  728. ! X(:,i) will be scaled to unit 2-norm
  729. WORK(i) = SCALE * ROOTSC
  730. CALL DLASCL( 'G',0, 0, WORK(i), ONE, M, 1, &
  731. Y(1,i), M, INFO2 ) ! LAPACK CALL
  732. ! Y(1:M,i) = (ONE/WORK(i)) * Y(1:M,i) ! INTRINSIC
  733. END IF
  734. ELSE
  735. WORK(i) = ZERO
  736. END IF
  737. END DO
  738. DO i = 1, N
  739. ! Now, apply the same scaling to the columns of X.
  740. IF ( WORK(i) > ZERO ) THEN
  741. CALL DSCAL( M, ONE/WORK(i), X(1,i), 1 ) ! BLAS CALL
  742. ! X(1:M,i) = (ONE/WORK(i)) * X(1:M,i) ! INTRINSIC
  743. ELSE IF ( WORK(i) < ZERO ) THEN
  744. CALL DLASCL( 'G', 0, 0, -WORK(i), &
  745. ONE/DBLE(M), M, 1, X(1,i), M, INFO2 ) ! LAPACK CALL
  746. ELSE IF ( X(IDAMAX(M, X(1,i),1),i ) &
  747. /= ZERO ) THEN
  748. ! Y(:,i) is zero vector. If X(:,i) is not
  749. ! zero, then a warning flag is raised.
  750. ! The computation continues but the
  751. ! situation will be reported in the output.
  752. BADXY = .TRUE.
  753. END IF
  754. END DO
  755. END IF
  756. !
  757. ! <2> SVD of the data snapshot matrix X.
  758. ! =====================================
  759. ! The left singular vectors are stored in the array X.
  760. ! The right singular vectors are in the array W.
  761. ! The array W will later on contain the eigenvectors
  762. ! of a Rayleigh quotient.
  763. NUMRNK = N
  764. SELECT CASE ( WHTSVD )
  765. CASE (1)
  766. CALL DGESVD( 'O', 'S', M, N, X, LDX, WORK, B, &
  767. LDB, W, LDW, WORK(N+1), LWORK-N, INFO1 ) ! LAPACK CALL
  768. T_OR_N = 'T'
  769. CASE (2)
  770. CALL DGESDD( 'O', M, N, X, LDX, WORK, B, LDB, W, &
  771. LDW, WORK(N+1), LWORK-N, IWORK, INFO1 ) ! LAPACK CALL
  772. T_OR_N = 'T'
  773. CASE (3)
  774. CALL DGESVDQ( 'H', 'P', 'N', 'R', 'R', M, N, &
  775. X, LDX, WORK, Z, LDZ, W, LDW, &
  776. NUMRNK, IWORK, LIWORK, WORK(N+MAX(2,M)+1),&
  777. LWORK-N-MAX(2,M), WORK(N+1), MAX(2,M), INFO1) ! LAPACK CALL
  778. CALL DLACPY( 'A', M, NUMRNK, Z, LDZ, X, LDX ) ! LAPACK CALL
  779. T_OR_N = 'T'
  780. CASE (4)
  781. CALL DGEJSV( 'F', 'U', JSVOPT, 'N', 'N', 'P', M, &
  782. N, X, LDX, WORK, Z, LDZ, W, LDW, &
  783. WORK(N+1), LWORK-N, IWORK, INFO1 ) ! LAPACK CALL
  784. CALL DLACPY( 'A', M, N, Z, LDZ, X, LDX ) ! LAPACK CALL
  785. T_OR_N = 'N'
  786. XSCL1 = WORK(N+1)
  787. XSCL2 = WORK(N+2)
  788. IF ( XSCL1 /= XSCL2 ) THEN
  789. ! This is an exceptional situation. If the
  790. ! data matrices are not scaled and the
  791. ! largest singular value of X overflows.
  792. ! In that case DGEJSV can return the SVD
  793. ! in scaled form. The scaling factor can be used
  794. ! to rescale the data (X and Y).
  795. CALL DLASCL( 'G', 0, 0, XSCL1, XSCL2, M, N, Y, LDY, INFO2 )
  796. END IF
  797. END SELECT
  798. !
  799. IF ( INFO1 > 0 ) THEN
  800. ! The SVD selected subroutine did not converge.
  801. ! Return with an error code.
  802. INFO = 2
  803. RETURN
  804. END IF
  805. !
  806. IF ( WORK(1) == ZERO ) THEN
  807. ! The largest computed singular value of (scaled)
  808. ! X is zero. Return error code -8
  809. ! (the 8th input variable had an illegal value).
  810. K = 0
  811. INFO = -8
  812. CALL XERBLA('DGEDMD',-INFO)
  813. RETURN
  814. END IF
  815. !
  816. !<3> Determine the numerical rank of the data
  817. ! snapshots matrix X. This depends on the
  818. ! parameters NRNK and TOL.
  819. SELECT CASE ( NRNK )
  820. CASE ( -1 )
  821. K = 1
  822. DO i = 2, NUMRNK
  823. IF ( ( WORK(i) <= WORK(1)*TOL ) .OR. &
  824. ( WORK(i) <= SMALL ) ) EXIT
  825. K = K + 1
  826. END DO
  827. CASE ( -2 )
  828. K = 1
  829. DO i = 1, NUMRNK-1
  830. IF ( ( WORK(i+1) <= WORK(i)*TOL ) .OR. &
  831. ( WORK(i) <= SMALL ) ) EXIT
  832. K = K + 1
  833. END DO
  834. CASE DEFAULT
  835. K = 1
  836. DO i = 2, NRNK
  837. IF ( WORK(i) <= SMALL ) EXIT
  838. K = K + 1
  839. END DO
  840. END SELECT
  841. ! Now, U = X(1:M,1:K) is the SVD/POD basis for the
  842. ! snapshot data in the input matrix X.
  843. !<4> Compute the Rayleigh quotient S = U^T * A * U.
  844. ! Depending on the requested outputs, the computation
  845. ! is organized to compute additional auxiliary
  846. ! matrices (for the residuals and refinements).
  847. !
  848. ! In all formulas below, we need V_k*Sigma_k^(-1)
  849. ! where either V_k is in W(1:N,1:K), or V_k^T is in
  850. ! W(1:K,1:N). Here Sigma_k=diag(WORK(1:K)).
  851. IF ( LSAME(T_OR_N, 'N') ) THEN
  852. DO i = 1, K
  853. CALL DSCAL( N, ONE/WORK(i), W(1,i), 1 ) ! BLAS CALL
  854. ! W(1:N,i) = (ONE/WORK(i)) * W(1:N,i) ! INTRINSIC
  855. END DO
  856. ELSE
  857. ! This non-unit stride access is due to the fact
  858. ! that DGESVD, DGESVDQ and DGESDD return the
  859. ! transposed matrix of the right singular vectors.
  860. !DO i = 1, K
  861. ! CALL DSCAL( N, ONE/WORK(i), W(i,1), LDW ) ! BLAS CALL
  862. ! ! W(i,1:N) = (ONE/WORK(i)) * W(i,1:N) ! INTRINSIC
  863. !END DO
  864. DO i = 1, K
  865. WORK(N+i) = ONE/WORK(i)
  866. END DO
  867. DO j = 1, N
  868. DO i = 1, K
  869. W(i,j) = (WORK(N+i))*W(i,j)
  870. END DO
  871. END DO
  872. END IF
  873. !
  874. IF ( WNTREF ) THEN
  875. !
  876. ! Need A*U(:,1:K)=Y*V_k*inv(diag(WORK(1:K)))
  877. ! for computing the refined Ritz vectors
  878. ! (optionally, outside DGEDMD).
  879. CALL DGEMM( 'N', T_OR_N, M, K, N, ONE, Y, LDY, W, &
  880. LDW, ZERO, Z, LDZ ) ! BLAS CALL
  881. ! Z(1:M,1:K)=MATMUL(Y(1:M,1:N),TRANSPOSE(W(1:K,1:N))) ! INTRINSIC, for T_OR_N=='T'
  882. ! Z(1:M,1:K)=MATMUL(Y(1:M,1:N),W(1:N,1:K)) ! INTRINSIC, for T_OR_N=='N'
  883. !
  884. ! At this point Z contains
  885. ! A * U(:,1:K) = Y * V_k * Sigma_k^(-1), and
  886. ! this is needed for computing the residuals.
  887. ! This matrix is returned in the array B and
  888. ! it can be used to compute refined Ritz vectors.
  889. CALL DLACPY( 'A', M, K, Z, LDZ, B, LDB ) ! BLAS CALL
  890. ! B(1:M,1:K) = Z(1:M,1:K) ! INTRINSIC
  891. CALL DGEMM( 'T', 'N', K, K, M, ONE, X, LDX, Z, &
  892. LDZ, ZERO, S, LDS ) ! BLAS CALL
  893. ! S(1:K,1:K) = MATMUL(TANSPOSE(X(1:M,1:K)),Z(1:M,1:K)) ! INTRINSIC
  894. ! At this point S = U^T * A * U is the Rayleigh quotient.
  895. ELSE
  896. ! A * U(:,1:K) is not explicitly needed and the
  897. ! computation is organized differently. The Rayleigh
  898. ! quotient is computed more efficiently.
  899. CALL DGEMM( 'T', 'N', K, N, M, ONE, X, LDX, Y, LDY, &
  900. ZERO, Z, LDZ ) ! BLAS CALL
  901. ! Z(1:K,1:N) = MATMUL( TRANSPOSE(X(1:M,1:K)), Y(1:M,1:N) ) ! INTRINSIC
  902. ! In the two DGEMM calls here, can use K for LDZ.
  903. CALL DGEMM( 'N', T_OR_N, K, K, N, ONE, Z, LDZ, W, &
  904. LDW, ZERO, S, LDS ) ! BLAS CALL
  905. ! S(1:K,1:K) = MATMUL(Z(1:K,1:N),TRANSPOSE(W(1:K,1:N))) ! INTRINSIC, for T_OR_N=='T'
  906. ! S(1:K,1:K) = MATMUL(Z(1:K,1:N),(W(1:N,1:K))) ! INTRINSIC, for T_OR_N=='N'
  907. ! At this point S = U^T * A * U is the Rayleigh quotient.
  908. ! If the residuals are requested, save scaled V_k into Z.
  909. ! Recall that V_k or V_k^T is stored in W.
  910. IF ( WNTRES .OR. WNTEX ) THEN
  911. IF ( LSAME(T_OR_N, 'N') ) THEN
  912. CALL DLACPY( 'A', N, K, W, LDW, Z, LDZ )
  913. ELSE
  914. CALL DLACPY( 'A', K, N, W, LDW, Z, LDZ )
  915. END IF
  916. END IF
  917. END IF
  918. !
  919. !<5> Compute the Ritz values and (if requested) the
  920. ! right eigenvectors of the Rayleigh quotient.
  921. !
  922. CALL DGEEV( 'N', JOBZL, K, S, LDS, REIG, IMEIG, W, &
  923. LDW, W, LDW, WORK(N+1), LWORK-N, INFO1 ) ! LAPACK CALL
  924. !
  925. ! W(1:K,1:K) contains the eigenvectors of the Rayleigh
  926. ! quotient. Even in the case of complex spectrum, all
  927. ! computation is done in real arithmetic. REIG and
  928. ! IMEIG are the real and the imaginary parts of the
  929. ! eigenvalues, so that the spectrum is given as
  930. ! REIG(:) + sqrt(-1)*IMEIG(:). Complex conjugate pairs
  931. ! are listed at consecutive positions. For such a
  932. ! complex conjugate pair of the eigenvalues, the
  933. ! corresponding eigenvectors are also a complex
  934. ! conjugate pair with the real and imaginary parts
  935. ! stored column-wise in W at the corresponding
  936. ! consecutive column indices. See the description of Z.
  937. ! Also, see the description of DGEEV.
  938. IF ( INFO1 > 0 ) THEN
  939. ! DGEEV failed to compute the eigenvalues and
  940. ! eigenvectors of the Rayleigh quotient.
  941. INFO = 3
  942. RETURN
  943. END IF
  944. !
  945. ! <6> Compute the eigenvectors (if requested) and,
  946. ! the residuals (if requested).
  947. !
  948. IF ( WNTVEC .OR. WNTEX ) THEN
  949. IF ( WNTRES ) THEN
  950. IF ( WNTREF ) THEN
  951. ! Here, if the refinement is requested, we have
  952. ! A*U(:,1:K) already computed and stored in Z.
  953. ! For the residuals, need Y = A * U(:,1;K) * W.
  954. CALL DGEMM( 'N', 'N', M, K, K, ONE, Z, LDZ, W, &
  955. LDW, ZERO, Y, LDY ) ! BLAS CALL
  956. ! Y(1:M,1:K) = Z(1:M,1:K) * W(1:K,1:K) ! INTRINSIC
  957. ! This frees Z; Y contains A * U(:,1:K) * W.
  958. ELSE
  959. ! Compute S = V_k * Sigma_k^(-1) * W, where
  960. ! V_k * Sigma_k^(-1) is stored in Z
  961. CALL DGEMM( T_OR_N, 'N', N, K, K, ONE, Z, LDZ, &
  962. W, LDW, ZERO, S, LDS)
  963. ! Then, compute Z = Y * S =
  964. ! = Y * V_k * Sigma_k^(-1) * W(1:K,1:K) =
  965. ! = A * U(:,1:K) * W(1:K,1:K)
  966. CALL DGEMM( 'N', 'N', M, K, N, ONE, Y, LDY, S, &
  967. LDS, ZERO, Z, LDZ)
  968. ! Save a copy of Z into Y and free Z for holding
  969. ! the Ritz vectors.
  970. CALL DLACPY( 'A', M, K, Z, LDZ, Y, LDY )
  971. IF ( WNTEX ) CALL DLACPY( 'A', M, K, Z, LDZ, B, LDB )
  972. END IF
  973. ELSE IF ( WNTEX ) THEN
  974. ! Compute S = V_k * Sigma_k^(-1) * W, where
  975. ! V_k * Sigma_k^(-1) is stored in Z
  976. CALL DGEMM( T_OR_N, 'N', N, K, K, ONE, Z, LDZ, &
  977. W, LDW, ZERO, S, LDS )
  978. ! Then, compute Z = Y * S =
  979. ! = Y * V_k * Sigma_k^(-1) * W(1:K,1:K) =
  980. ! = A * U(:,1:K) * W(1:K,1:K)
  981. CALL DGEMM( 'N', 'N', M, K, N, ONE, Y, LDY, S, &
  982. LDS, ZERO, B, LDB )
  983. ! The above call replaces the following two calls
  984. ! that were used in the developing-testing phase.
  985. ! CALL DGEMM( 'N', 'N', M, K, N, ONE, Y, LDY, S, &
  986. ! LDS, ZERO, Z, LDZ)
  987. ! Save a copy of Z into B and free Z for holding
  988. ! the Ritz vectors.
  989. ! CALL DLACPY( 'A', M, K, Z, LDZ, B, LDB )
  990. END IF
  991. !
  992. ! Compute the real form of the Ritz vectors
  993. IF ( WNTVEC ) CALL DGEMM( 'N', 'N', M, K, K, ONE, X, LDX, W, LDW, &
  994. ZERO, Z, LDZ ) ! BLAS CALL
  995. ! Z(1:M,1:K) = MATMUL(X(1:M,1:K), W(1:K,1:K)) ! INTRINSIC
  996. !
  997. IF ( WNTRES ) THEN
  998. i = 1
  999. DO WHILE ( i <= K )
  1000. IF ( IMEIG(i) == ZERO ) THEN
  1001. ! have a real eigenvalue with real eigenvector
  1002. CALL DAXPY( M, -REIG(i), Z(1,i), 1, Y(1,i), 1 ) ! BLAS CALL
  1003. ! Y(1:M,i) = Y(1:M,i) - REIG(i) * Z(1:M,i) ! INTRINSIC
  1004. RES(i) = DNRM2( M, Y(1,i), 1) ! BLAS CALL
  1005. i = i + 1
  1006. ELSE
  1007. ! Have a complex conjugate pair
  1008. ! REIG(i) +- sqrt(-1)*IMEIG(i).
  1009. ! Since all computation is done in real
  1010. ! arithmetic, the formula for the residual
  1011. ! is recast for real representation of the
  1012. ! complex conjugate eigenpair. See the
  1013. ! description of RES.
  1014. AB(1,1) = REIG(i)
  1015. AB(2,1) = -IMEIG(i)
  1016. AB(1,2) = IMEIG(i)
  1017. AB(2,2) = REIG(i)
  1018. CALL DGEMM( 'N', 'N', M, 2, 2, -ONE, Z(1,i), &
  1019. LDZ, AB, 2, ONE, Y(1,i), LDY ) ! BLAS CALL
  1020. ! Y(1:M,i:i+1) = Y(1:M,i:i+1) - Z(1:M,i:i+1) * AB ! INTRINSIC
  1021. RES(i) = DLANGE( 'F', M, 2, Y(1,i), LDY, &
  1022. WORK(N+1) ) ! LAPACK CALL
  1023. RES(i+1) = RES(i)
  1024. i = i + 2
  1025. END IF
  1026. END DO
  1027. END IF
  1028. END IF
  1029. !
  1030. IF ( WHTSVD == 4 ) THEN
  1031. WORK(N+1) = XSCL1
  1032. WORK(N+2) = XSCL2
  1033. END IF
  1034. !
  1035. ! Successful exit.
  1036. IF ( .NOT. BADXY ) THEN
  1037. INFO = 0
  1038. ELSE
  1039. ! A warning on possible data inconsistency.
  1040. ! This should be a rare event.
  1041. INFO = 4
  1042. END IF
  1043. !............................................................
  1044. RETURN
  1045. ! ......
  1046. END SUBROUTINE DGEDMD