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.

sgedmd.f90 46 kB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054
  1. SUBROUTINE SGEDMD( 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 = real32
  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. ! SGEDMD 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, SGEDMD 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, SGEDMD 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 :: SGESVD (the QR SVD algorithm)
  140. ! 2 :: SGESDD (the Divide and Conquer algorithm; if enough
  141. ! workspace available, this is the fastest option)
  142. ! 3 :: SGESVDQ (the preconditioned QR SVD ; this and 4
  143. ! are the most accurate options)
  144. ! 4 :: SGEJSV (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. ! The numerical rank can be enforced by using positive
  202. ! value of NRNK as follows:
  203. ! 0 < NRNK <= N :: at most NRNK largest singular values
  204. ! will be used. If the number of the computed nonzero
  205. ! singular values is less than NRNK, then only those
  206. ! nonzero values will be used and the actually used
  207. ! dimension is less than NRNK. The actual number of
  208. ! the nonzero singular values is returned in the variable
  209. ! K. See the descriptions of TOL and K.
  210. !.....
  211. ! TOL (input) REAL(KIND=WP), 0 <= TOL < 1
  212. ! The tolerance for truncating small singular values.
  213. ! See the description of NRNK.
  214. !.....
  215. ! K (output) INTEGER, 0 <= K <= N
  216. ! The dimension of the POD basis for the data snapshot
  217. ! matrix X and the number of the computed Ritz pairs.
  218. ! The value of K is determined according to the rule set
  219. ! by the parameters NRNK and TOL.
  220. ! See the descriptions of NRNK and TOL.
  221. !.....
  222. ! REIG (output) REAL(KIND=WP) N-by-1 array
  223. ! The leading K (K<=N) entries of REIG contain
  224. ! the real parts of the computed eigenvalues
  225. ! REIG(1:K) + sqrt(-1)*IMEIG(1:K).
  226. ! See the descriptions of K, IMEIG, and Z.
  227. !.....
  228. ! IMEIG (output) REAL(KIND=WP) N-by-1 array
  229. ! The leading K (K<=N) entries of IMEIG contain
  230. ! the imaginary parts of the computed eigenvalues
  231. ! REIG(1:K) + sqrt(-1)*IMEIG(1:K).
  232. ! The eigenvalues are determined as follows:
  233. ! If IMEIG(i) == 0, then the corresponding eigenvalue is
  234. ! real, LAMBDA(i) = REIG(i).
  235. ! If IMEIG(i)>0, then the corresponding complex
  236. ! conjugate pair of eigenvalues reads
  237. ! LAMBDA(i) = REIG(i) + sqrt(-1)*IMAG(i)
  238. ! LAMBDA(i+1) = REIG(i) - sqrt(-1)*IMAG(i)
  239. ! That is, complex conjugate pairs have consecutive
  240. ! indices (i,i+1), with the positive imaginary part
  241. ! listed first.
  242. ! See the descriptions of K, REIG, and Z.
  243. !.....
  244. ! Z (workspace/output) REAL(KIND=WP) M-by-N array
  245. ! If JOBZ =='V' then
  246. ! Z contains real Ritz vectors as follows:
  247. ! If IMEIG(i)=0, then Z(:,i) is an eigenvector of
  248. ! the i-th Ritz value; ||Z(:,i)||_2=1.
  249. ! If IMEIG(i) > 0 (and IMEIG(i+1) < 0) then
  250. ! [Z(:,i) Z(:,i+1)] span an invariant subspace and
  251. ! the Ritz values extracted from this subspace are
  252. ! REIG(i) + sqrt(-1)*IMEIG(i) and
  253. ! REIG(i) - sqrt(-1)*IMEIG(i).
  254. ! The corresponding eigenvectors are
  255. ! Z(:,i) + sqrt(-1)*Z(:,i+1) and
  256. ! Z(:,i) - sqrt(-1)*Z(:,i+1), respectively.
  257. ! || Z(:,i:i+1)||_F = 1.
  258. ! If JOBZ == 'F', then the above descriptions hold for
  259. ! the columns of X(:,1:K)*W(1:K,1:K), where the columns
  260. ! of W(1:k,1:K) are the computed eigenvectors of the
  261. ! K-by-K Rayleigh quotient. The columns of W(1:K,1:K)
  262. ! are similarly structured: If IMEIG(i) == 0 then
  263. ! X(:,1:K)*W(:,i) is an eigenvector, and if IMEIG(i)>0
  264. ! then X(:,1:K)*W(:,i)+sqrt(-1)*X(:,1:K)*W(:,i+1) and
  265. ! X(:,1:K)*W(:,i)-sqrt(-1)*X(:,1:K)*W(:,i+1)
  266. ! are the eigenvectors of LAMBDA(i), LAMBDA(i+1).
  267. ! See the descriptions of REIG, IMEIG, X and W.
  268. !.....
  269. ! LDZ (input) INTEGER , LDZ >= M
  270. ! The leading dimension of the array Z.
  271. !.....
  272. ! RES (output) REAL(KIND=WP) N-by-1 array
  273. ! RES(1:K) contains the residuals for the K computed
  274. ! Ritz pairs.
  275. ! If LAMBDA(i) is real, then
  276. ! RES(i) = || A * Z(:,i) - LAMBDA(i)*Z(:,i))||_2.
  277. ! If [LAMBDA(i), LAMBDA(i+1)] is a complex conjugate pair
  278. ! then
  279. ! RES(i)=RES(i+1) = || A * Z(:,i:i+1) - Z(:,i:i+1) *B||_F
  280. ! where B = [ real(LAMBDA(i)) imag(LAMBDA(i)) ]
  281. ! [-imag(LAMBDA(i)) real(LAMBDA(i)) ].
  282. ! It holds that
  283. ! RES(i) = || A*ZC(:,i) - LAMBDA(i) *ZC(:,i) ||_2
  284. ! RES(i+1) = || A*ZC(:,i+1) - LAMBDA(i+1)*ZC(:,i+1) ||_2
  285. ! where ZC(:,i) = Z(:,i) + sqrt(-1)*Z(:,i+1)
  286. ! ZC(:,i+1) = Z(:,i) - sqrt(-1)*Z(:,i+1)
  287. ! See the description of REIG, IMEIG and Z.
  288. !.....
  289. ! B (output) REAL(KIND=WP) M-by-N array.
  290. ! IF JOBF =='R', B(1:M,1:K) contains A*U(:,1:K), and can
  291. ! be used for computing the refined vectors; see further
  292. ! details in the provided references.
  293. ! If JOBF == 'E', B(1:M,1;K) contains
  294. ! A*U(:,1:K)*W(1:K,1:K), which are the vectors from the
  295. ! Exact DMD, up to scaling by the inverse eigenvalues.
  296. ! If JOBF =='N', then B is not referenced.
  297. ! See the descriptions of X, W, K.
  298. !.....
  299. ! LDB (input) INTEGER, LDB >= M
  300. ! The leading dimension of the array B.
  301. !.....
  302. ! W (workspace/output) REAL(KIND=WP) N-by-N array
  303. ! On exit, W(1:K,1:K) contains the K computed
  304. ! eigenvectors of the matrix Rayleigh quotient (real and
  305. ! imaginary parts for each complex conjugate pair of the
  306. ! eigenvalues). The Ritz vectors (returned in Z) are the
  307. ! product of X (containing a POD basis for the input
  308. ! matrix X) and W. See the descriptions of K, S, X and Z.
  309. ! W is also used as a workspace to temporarily store the
  310. ! left singular vectors of X.
  311. !.....
  312. ! LDW (input) INTEGER, LDW >= N
  313. ! The leading dimension of the array W.
  314. !.....
  315. ! S (workspace/output) REAL(KIND=WP) N-by-N array
  316. ! The array S(1:K,1:K) is used for the matrix Rayleigh
  317. ! quotient. This content is overwritten during
  318. ! the eigenvalue decomposition by SGEEV.
  319. ! See the description of K.
  320. !.....
  321. ! LDS (input) INTEGER, LDS >= N
  322. ! The leading dimension of the array S.
  323. !.....
  324. ! WORK (workspace/output) REAL(KIND=WP) LWORK-by-1 array
  325. ! On exit, WORK(1:N) contains the singular values of
  326. ! X (for JOBS=='N') or column scaled X (JOBS=='S', 'C').
  327. ! If WHTSVD==4, then WORK(N+1) and WORK(N+2) contain
  328. ! scaling factor WORK(N+2)/WORK(N+1) used to scale X
  329. ! and Y to avoid overflow in the SVD of X.
  330. ! This may be of interest if the scaling option is off
  331. ! and as many as possible smallest eigenvalues are
  332. ! desired to the highest feasible accuracy.
  333. ! If the call to SGEDMD is only workspace query, then
  334. ! WORK(1) contains the minimal workspace length and
  335. ! WORK(2) is the optimal workspace length. Hence, the
  336. ! length of work is at least 2.
  337. ! See the description of LWORK.
  338. !.....
  339. ! LWORK (input) INTEGER
  340. ! The minimal length of the workspace vector WORK.
  341. ! LWORK is calculated as follows:
  342. ! If WHTSVD == 1 ::
  343. ! If JOBZ == 'V', then
  344. ! LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,4*N)).
  345. ! If JOBZ == 'N' then
  346. ! LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,3*N)).
  347. ! Here LWORK_SVD = MAX(1,3*N+M,5*N) is the minimal
  348. ! workspace length of SGESVD.
  349. ! If WHTSVD == 2 ::
  350. ! If JOBZ == 'V', then
  351. ! LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,4*N))
  352. ! If JOBZ == 'N', then
  353. ! LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,3*N))
  354. ! Here LWORK_SVD = MAX(M, 5*N*N+4*N)+3*N*N is the
  355. ! minimal workspace length of SGESDD.
  356. ! If WHTSVD == 3 ::
  357. ! If JOBZ == 'V', then
  358. ! LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,4*N))
  359. ! If JOBZ == 'N', then
  360. ! LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,3*N))
  361. ! Here LWORK_SVD = N+M+MAX(3*N+1,
  362. ! MAX(1,3*N+M,5*N),MAX(1,N))
  363. ! is the minimal workspace length of SGESVDQ.
  364. ! If WHTSVD == 4 ::
  365. ! If JOBZ == 'V', then
  366. ! LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,4*N))
  367. ! If JOBZ == 'N', then
  368. ! LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,3*N))
  369. ! Here LWORK_SVD = MAX(7,2*M+N,6*N+2*N*N) is the
  370. ! minimal workspace length of SGEJSV.
  371. ! The above expressions are not simplified in order to
  372. ! make the usage of WORK more transparent, and for
  373. ! easier checking. In any case, LWORK >= 2.
  374. ! If on entry LWORK = -1, then a workspace query is
  375. ! assumed and the procedure only computes the minimal
  376. ! and the optimal workspace lengths for both WORK and
  377. ! IWORK. See the descriptions of WORK and IWORK.
  378. !.....
  379. ! IWORK (workspace/output) INTEGER LIWORK-by-1 array
  380. ! Workspace that is required only if WHTSVD equals
  381. ! 2 , 3 or 4. (See the description of WHTSVD).
  382. ! If on entry LWORK =-1 or LIWORK=-1, then the
  383. ! minimal length of IWORK is computed and returned in
  384. ! IWORK(1). See the description of LIWORK.
  385. !.....
  386. ! LIWORK (input) INTEGER
  387. ! The minimal length of the workspace vector IWORK.
  388. ! If WHTSVD == 1, then only IWORK(1) is used; LIWORK >=1
  389. ! If WHTSVD == 2, then LIWORK >= MAX(1,8*MIN(M,N))
  390. ! If WHTSVD == 3, then LIWORK >= MAX(1,M+N-1)
  391. ! If WHTSVD == 4, then LIWORK >= MAX(3,M+3*N)
  392. ! If on entry LIWORK = -1, then a workspace query is
  393. ! assumed and the procedure only computes the minimal
  394. ! and the optimal workspace lengths for both WORK and
  395. ! IWORK. See the descriptions of WORK and IWORK.
  396. !.....
  397. ! INFO (output) INTEGER
  398. ! -i < 0 :: On entry, the i-th argument had an
  399. ! illegal value
  400. ! = 0 :: Successful return.
  401. ! = 1 :: Void input. Quick exit (M=0 or N=0).
  402. ! = 2 :: The SVD computation of X did not converge.
  403. ! Suggestion: Check the input data and/or
  404. ! repeat with different WHTSVD.
  405. ! = 3 :: The computation of the eigenvalues did not
  406. ! converge.
  407. ! = 4 :: If data scaling was requested on input and
  408. ! the procedure found inconsistency in the data
  409. ! such that for some column index i,
  410. ! X(:,i) = 0 but Y(:,i) /= 0, then Y(:,i) is set
  411. ! to zero if JOBS=='C'. The computation proceeds
  412. ! with original or modified data and warning
  413. ! flag is set with INFO=4.
  414. !.............................................................
  415. !.............................................................
  416. ! Parameters
  417. ! ~~~~~~~~~~
  418. REAL(KIND=WP), PARAMETER :: ONE = 1.0_WP
  419. REAL(KIND=WP), PARAMETER :: ZERO = 0.0_WP
  420. ! Local scalars
  421. ! ~~~~~~~~~~~~~
  422. REAL(KIND=WP) :: OFL, ROOTSC, SCALE, SMALL, &
  423. SSUM, XSCL1, XSCL2
  424. INTEGER :: i, j, IMINWR, INFO1, INFO2, &
  425. LWRKEV, LWRSDD, LWRSVD, &
  426. LWRSVQ, MLWORK, MWRKEV, MWRSDD, &
  427. MWRSVD, MWRSVJ, MWRSVQ, NUMRNK, &
  428. OLWORK
  429. LOGICAL :: BADXY, LQUERY, SCCOLX, SCCOLY, &
  430. WNTEX, WNTREF, WNTRES, WNTVEC
  431. CHARACTER :: JOBZL, T_OR_N
  432. CHARACTER :: JSVOPT
  433. ! Local arrays
  434. ! ~~~~~~~~~~~~
  435. REAL(KIND=WP) :: AB(2,2), RDUMMY(2), RDUMMY2(2)
  436. ! External functions (BLAS and LAPACK)
  437. ! ~~~~~~~~~~~~~~~~~
  438. REAL(KIND=WP) SLANGE, SLAMCH, SNRM2
  439. EXTERNAL SLANGE, SLAMCH, SNRM2, ISAMAX
  440. INTEGER ISAMAX
  441. LOGICAL SISNAN, LSAME
  442. EXTERNAL SISNAN, LSAME
  443. ! External subroutines (BLAS and LAPACK)
  444. ! ~~~~~~~~~~~~~~~~~~~~
  445. EXTERNAL SAXPY, SGEMM, SSCAL
  446. EXTERNAL SGEEV, SGEJSV, SGESDD, SGESVD, SGESVDQ, &
  447. SLACPY, SLASCL, SLASSQ, XERBLA
  448. ! Intrinsic functions
  449. ! ~~~~~~~~~~~~~~~~~~~
  450. INTRINSIC INT, FLOAT, MAX, SQRT
  451. !............................................................
  452. !
  453. ! Test the input arguments
  454. !
  455. WNTRES = LSAME(JOBR,'R')
  456. SCCOLX = LSAME(JOBS,'S') .OR. LSAME(JOBS,'C')
  457. SCCOLY = LSAME(JOBS,'Y')
  458. WNTVEC = LSAME(JOBZ,'V')
  459. WNTREF = LSAME(JOBF,'R')
  460. WNTEX = LSAME(JOBF,'E')
  461. INFO = 0
  462. LQUERY = ( ( LWORK == -1 ) .OR. ( LIWORK == -1 ) )
  463. !
  464. IF ( .NOT. (SCCOLX .OR. SCCOLY .OR. &
  465. LSAME(JOBS,'N')) ) THEN
  466. INFO = -1
  467. ELSE IF ( .NOT. (WNTVEC .OR. LSAME(JOBZ,'N') &
  468. .OR. LSAME(JOBZ,'F')) ) THEN
  469. INFO = -2
  470. ELSE IF ( .NOT. (WNTRES .OR. LSAME(JOBR,'N')) .OR. &
  471. ( WNTRES .AND. (.NOT.WNTVEC) ) ) THEN
  472. INFO = -3
  473. ELSE IF ( .NOT. (WNTREF .OR. WNTEX .OR. &
  474. LSAME(JOBF,'N') ) ) THEN
  475. INFO = -4
  476. ELSE IF ( .NOT.((WHTSVD == 1) .OR. (WHTSVD == 2) .OR. &
  477. (WHTSVD == 3) .OR. (WHTSVD == 4) )) THEN
  478. INFO = -5
  479. ELSE IF ( M < 0 ) THEN
  480. INFO = -6
  481. ELSE IF ( ( N < 0 ) .OR. ( N > M ) ) THEN
  482. INFO = -7
  483. ELSE IF ( LDX < M ) THEN
  484. INFO = -9
  485. ELSE IF ( LDY < M ) THEN
  486. INFO = -11
  487. ELSE IF ( .NOT. (( NRNK == -2).OR.(NRNK == -1).OR. &
  488. ((NRNK >= 1).AND.(NRNK <=N ))) ) THEN
  489. INFO = -12
  490. ELSE IF ( ( TOL < ZERO ) .OR. ( TOL >= ONE ) ) THEN
  491. INFO = -13
  492. ELSE IF ( LDZ < M ) THEN
  493. INFO = -18
  494. ELSE IF ( (WNTREF .OR. WNTEX ) .AND. ( LDB < M ) ) THEN
  495. INFO = -21
  496. ELSE IF ( LDW < N ) THEN
  497. INFO = -23
  498. ELSE IF ( LDS < N ) THEN
  499. INFO = -25
  500. END IF
  501. !
  502. IF ( INFO == 0 ) THEN
  503. ! Compute the minimal and the optimal workspace
  504. ! requirements. Simulate running the code and
  505. ! determine minimal and optimal sizes of the
  506. ! workspace at any moment of the run.
  507. IF ( N == 0 ) THEN
  508. ! Quick return. All output except K is void.
  509. ! INFO=1 signals the void input.
  510. ! In case of a workspace query, the default
  511. ! minimal workspace lengths are returned.
  512. IF ( LQUERY ) THEN
  513. IWORK(1) = 1
  514. WORK(1) = 2
  515. WORK(2) = 2
  516. ELSE
  517. K = 0
  518. END IF
  519. INFO = 1
  520. RETURN
  521. END IF
  522. MLWORK = MAX(2,N)
  523. OLWORK = MAX(2,N)
  524. IMINWR = 1
  525. SELECT CASE ( WHTSVD )
  526. CASE (1)
  527. ! The following is specified as the minimal
  528. ! length of WORK in the definition of SGESVD:
  529. ! MWRSVD = MAX(1,3*MIN(M,N)+MAX(M,N),5*MIN(M,N))
  530. MWRSVD = MAX(1,3*MIN(M,N)+MAX(M,N),5*MIN(M,N))
  531. MLWORK = MAX(MLWORK,N + MWRSVD)
  532. IF ( LQUERY ) THEN
  533. CALL SGESVD( 'O', 'S', M, N, X, LDX, WORK, &
  534. B, LDB, W, LDW, RDUMMY, -1, INFO1 )
  535. LWRSVD = MAX( MWRSVD, INT( RDUMMY(1) ) )
  536. OLWORK = MAX(OLWORK,N + LWRSVD)
  537. END IF
  538. CASE (2)
  539. ! The following is specified as the minimal
  540. ! length of WORK in the definition of SGESDD:
  541. ! MWRSDD = 3*MIN(M,N)*MIN(M,N) +
  542. ! MAX( MAX(M,N),5*MIN(M,N)*MIN(M,N)+4*MIN(M,N) )
  543. ! IMINWR = 8*MIN(M,N)
  544. MWRSDD = 3*MIN(M,N)*MIN(M,N) + &
  545. MAX( MAX(M,N),5*MIN(M,N)*MIN(M,N)+4*MIN(M,N) )
  546. MLWORK = MAX(MLWORK,N + MWRSDD)
  547. IMINWR = 8*MIN(M,N)
  548. IF ( LQUERY ) THEN
  549. CALL SGESDD( 'O', M, N, X, LDX, WORK, B, &
  550. LDB, W, LDW, RDUMMY, -1, IWORK, INFO1 )
  551. LWRSDD = MAX( MWRSDD, INT( RDUMMY(1) ) )
  552. OLWORK = MAX(OLWORK,N + LWRSDD)
  553. END IF
  554. CASE (3)
  555. !LWQP3 = 3*N+1
  556. !LWORQ = MAX(N, 1)
  557. !MWRSVD = MAX(1,3*MIN(M,N)+MAX(M,N),5*MIN(M,N))
  558. !MWRSVQ = N + MAX( LWQP3, MWRSVD, LWORQ )+ MAX(M,2)
  559. !MLWORK = N + MWRSVQ
  560. !IMINWR = M+N-1
  561. CALL SGESVDQ( 'H', 'P', 'N', 'R', 'R', M, N, &
  562. X, LDX, WORK, Z, LDZ, W, LDW, &
  563. NUMRNK, IWORK, -1, RDUMMY, &
  564. -1, RDUMMY2, -1, INFO1 )
  565. IMINWR = IWORK(1)
  566. MWRSVQ = INT(RDUMMY(2))
  567. MLWORK = MAX(MLWORK,N+MWRSVQ+INT(RDUMMY2(1)))
  568. IF ( LQUERY ) THEN
  569. LWRSVQ = INT(RDUMMY(1))
  570. OLWORK = MAX(OLWORK,N+LWRSVQ+INT(RDUMMY2(1)))
  571. END IF
  572. CASE (4)
  573. JSVOPT = 'J'
  574. !MWRSVJ = MAX( 7, 2*M+N, 6*N+2*N*N )! for JSVOPT='V'
  575. MWRSVJ = MAX( 7, 2*M+N, 4*N+N*N, 2*N+N*N+6 )
  576. MLWORK = MAX(MLWORK,N+MWRSVJ)
  577. IMINWR = MAX( 3, M+3*N )
  578. IF ( LQUERY ) THEN
  579. OLWORK = MAX(OLWORK,N+MWRSVJ)
  580. END IF
  581. END SELECT
  582. IF ( WNTVEC .OR. WNTEX .OR. LSAME(JOBZ,'F') ) THEN
  583. JOBZL = 'V'
  584. ELSE
  585. JOBZL = 'N'
  586. END IF
  587. ! Workspace calculation to the SGEEV call
  588. IF ( LSAME(JOBZL,'V') ) THEN
  589. MWRKEV = MAX( 1, 4*N )
  590. ELSE
  591. MWRKEV = MAX( 1, 3*N )
  592. END IF
  593. MLWORK = MAX(MLWORK,N+MWRKEV)
  594. IF ( LQUERY ) THEN
  595. CALL SGEEV( 'N', JOBZL, N, S, LDS, REIG, &
  596. IMEIG, W, LDW, W, LDW, RDUMMY, -1, INFO1 )
  597. LWRKEV = MAX( MWRKEV, INT(RDUMMY(1)) )
  598. OLWORK = MAX( OLWORK, N+LWRKEV )
  599. END IF
  600. !
  601. IF ( LIWORK < IMINWR .AND. (.NOT.LQUERY) ) INFO = -29
  602. IF ( LWORK < MLWORK .AND. (.NOT.LQUERY) ) INFO = -27
  603. END IF
  604. !
  605. IF( INFO /= 0 ) THEN
  606. CALL XERBLA( 'SGEDMD', -INFO )
  607. RETURN
  608. ELSE IF ( LQUERY ) THEN
  609. ! Return minimal and optimal workspace sizes
  610. IWORK(1) = IMINWR
  611. WORK(1) = MLWORK
  612. WORK(2) = OLWORK
  613. RETURN
  614. END IF
  615. !............................................................
  616. !
  617. OFL = SLAMCH('O')
  618. SMALL = SLAMCH('S')
  619. BADXY = .FALSE.
  620. !
  621. ! <1> Optional scaling of the snapshots (columns of X, Y)
  622. ! ==========================================================
  623. IF ( SCCOLX ) THEN
  624. ! The columns of X will be normalized.
  625. ! To prevent overflows, the column norms of X are
  626. ! carefully computed using SLASSQ.
  627. K = 0
  628. DO i = 1, N
  629. !WORK(i) = DNRM2( M, X(1,i), 1 )
  630. SCALE = ZERO
  631. CALL SLASSQ( M, X(1,i), 1, SCALE, SSUM )
  632. IF ( SISNAN(SCALE) .OR. SISNAN(SSUM) ) THEN
  633. K = 0
  634. INFO = -8
  635. CALL XERBLA('SGEDMD',-INFO)
  636. END IF
  637. IF ( (SCALE /= ZERO) .AND. (SSUM /= ZERO) ) THEN
  638. ROOTSC = SQRT(SSUM)
  639. IF ( SCALE .GE. (OFL / ROOTSC) ) THEN
  640. ! Norm of X(:,i) overflows. First, X(:,i)
  641. ! is scaled by
  642. ! ( ONE / ROOTSC ) / SCALE = 1/||X(:,i)||_2.
  643. ! Next, the norm of X(:,i) is stored without
  644. ! overflow as WORK(i) = - SCALE * (ROOTSC/M),
  645. ! the minus sign indicating the 1/M factor.
  646. ! Scaling is performed without overflow, and
  647. ! underflow may occur in the smallest entries
  648. ! of X(:,i). The relative backward and forward
  649. ! errors are small in the ell_2 norm.
  650. CALL SLASCL( 'G', 0, 0, SCALE, ONE/ROOTSC, &
  651. M, 1, X(1,i), M, INFO2 )
  652. WORK(i) = - SCALE * ( ROOTSC / FLOAT(M) )
  653. ELSE
  654. ! X(:,i) will be scaled to unit 2-norm
  655. WORK(i) = SCALE * ROOTSC
  656. CALL SLASCL( 'G',0, 0, WORK(i), ONE, M, 1, &
  657. X(1,i), M, INFO2 ) ! LAPACK CALL
  658. ! X(1:M,i) = (ONE/WORK(i)) * X(1:M,i) ! INTRINSIC
  659. END IF
  660. ELSE
  661. WORK(i) = ZERO
  662. K = K + 1
  663. END IF
  664. END DO
  665. IF ( K == N ) THEN
  666. ! All columns of X are zero. Return error code -8.
  667. ! (the 8th input variable had an illegal value)
  668. K = 0
  669. INFO = -8
  670. CALL XERBLA('SGEDMD',-INFO)
  671. RETURN
  672. END IF
  673. DO i = 1, N
  674. ! Now, apply the same scaling to the columns of Y.
  675. IF ( WORK(i) > ZERO ) THEN
  676. CALL SSCAL( M, ONE/WORK(i), Y(1,i), 1 ) ! BLAS CALL
  677. ! Y(1:M,i) = (ONE/WORK(i)) * Y(1:M,i) ! INTRINSIC
  678. ELSE IF ( WORK(i) < ZERO ) THEN
  679. CALL SLASCL( 'G', 0, 0, -WORK(i), &
  680. ONE/FLOAT(M), M, 1, Y(1,i), M, INFO2 ) ! LAPACK CALL
  681. ELSE IF ( Y(ISAMAX(M, Y(1,i),1),i ) &
  682. /= ZERO ) THEN
  683. ! X(:,i) is zero vector. For consistency,
  684. ! Y(:,i) should also be zero. If Y(:,i) is not
  685. ! zero, then the data might be inconsistent or
  686. ! corrupted. If JOBS == 'C', Y(:,i) is set to
  687. ! zero and a warning flag is raised.
  688. ! The computation continues but the
  689. ! situation will be reported in the output.
  690. BADXY = .TRUE.
  691. IF ( LSAME(JOBS,'C')) &
  692. CALL SSCAL( M, ZERO, Y(1,i), 1 ) ! BLAS CALL
  693. END IF
  694. END DO
  695. END IF
  696. !
  697. IF ( SCCOLY ) THEN
  698. ! The columns of Y will be normalized.
  699. ! To prevent overflows, the column norms of Y are
  700. ! carefully computed using SLASSQ.
  701. DO i = 1, N
  702. !WORK(i) = DNRM2( M, Y(1,i), 1 )
  703. SCALE = ZERO
  704. CALL SLASSQ( M, Y(1,i), 1, SCALE, SSUM )
  705. IF ( SISNAN(SCALE) .OR. SISNAN(SSUM) ) THEN
  706. K = 0
  707. INFO = -10
  708. CALL XERBLA('SGEDMD',-INFO)
  709. END IF
  710. IF ( SCALE /= ZERO .AND. (SSUM /= ZERO) ) THEN
  711. ROOTSC = SQRT(SSUM)
  712. IF ( SCALE .GE. (OFL / ROOTSC) ) THEN
  713. ! Norm of Y(:,i) overflows. First, Y(:,i)
  714. ! is scaled by
  715. ! ( ONE / ROOTSC ) / SCALE = 1/||Y(:,i)||_2.
  716. ! Next, the norm of Y(:,i) is stored without
  717. ! overflow as WORK(i) = - SCALE * (ROOTSC/M),
  718. ! the minus sign indicating the 1/M factor.
  719. ! Scaling is performed without overflow, and
  720. ! underflow may occur in the smallest entries
  721. ! of Y(:,i). The relative backward and forward
  722. ! errors are small in the ell_2 norm.
  723. CALL SLASCL( 'G', 0, 0, SCALE, ONE/ROOTSC, &
  724. M, 1, Y(1,i), M, INFO2 )
  725. WORK(i) = - SCALE * ( ROOTSC / FLOAT(M) )
  726. ELSE
  727. ! X(:,i) will be scaled to unit 2-norm
  728. WORK(i) = SCALE * ROOTSC
  729. CALL SLASCL( 'G',0, 0, WORK(i), ONE, M, 1, &
  730. Y(1,i), M, INFO2 ) ! LAPACK CALL
  731. ! Y(1:M,i) = (ONE/WORK(i)) * Y(1:M,i) ! INTRINSIC
  732. END IF
  733. ELSE
  734. WORK(i) = ZERO
  735. END IF
  736. END DO
  737. DO i = 1, N
  738. ! Now, apply the same scaling to the columns of X.
  739. IF ( WORK(i) > ZERO ) THEN
  740. CALL SSCAL( M, ONE/WORK(i), X(1,i), 1 ) ! BLAS CALL
  741. ! X(1:M,i) = (ONE/WORK(i)) * X(1:M,i) ! INTRINSIC
  742. ELSE IF ( WORK(i) < ZERO ) THEN
  743. CALL SLASCL( 'G', 0, 0, -WORK(i), &
  744. ONE/FLOAT(M), M, 1, X(1,i), M, INFO2 ) ! LAPACK CALL
  745. ELSE IF ( X(ISAMAX(M, X(1,i),1),i ) &
  746. /= ZERO ) THEN
  747. ! Y(:,i) is zero vector. If X(:,i) is not
  748. ! zero, then a warning flag is raised.
  749. ! The computation continues but the
  750. ! situation will be reported in the output.
  751. BADXY = .TRUE.
  752. END IF
  753. END DO
  754. END IF
  755. !
  756. ! <2> SVD of the data snapshot matrix X.
  757. ! =====================================
  758. ! The left singular vectors are stored in the array X.
  759. ! The right singular vectors are in the array W.
  760. ! The array W will later on contain the eigenvectors
  761. ! of a Rayleigh quotient.
  762. NUMRNK = N
  763. SELECT CASE ( WHTSVD )
  764. CASE (1)
  765. CALL SGESVD( 'O', 'S', M, N, X, LDX, WORK, B, &
  766. LDB, W, LDW, WORK(N+1), LWORK-N, INFO1 ) ! LAPACK CALL
  767. T_OR_N = 'T'
  768. CASE (2)
  769. CALL SGESDD( 'O', M, N, X, LDX, WORK, B, LDB, W, &
  770. LDW, WORK(N+1), LWORK-N, IWORK, INFO1 ) ! LAPACK CALL
  771. T_OR_N = 'T'
  772. CASE (3)
  773. CALL SGESVDQ( 'H', 'P', 'N', 'R', 'R', M, N, &
  774. X, LDX, WORK, Z, LDZ, W, LDW, &
  775. NUMRNK, IWORK, LIWORK, WORK(N+MAX(2,M)+1),&
  776. LWORK-N-MAX(2,M), WORK(N+1), MAX(2,M), INFO1) ! LAPACK CALL
  777. CALL SLACPY( 'A', M, NUMRNK, Z, LDZ, X, LDX ) ! LAPACK CALL
  778. T_OR_N = 'T'
  779. CASE (4)
  780. CALL SGEJSV( 'F', 'U', JSVOPT, 'N', 'N', 'P', M, &
  781. N, X, LDX, WORK, Z, LDZ, W, LDW, &
  782. WORK(N+1), LWORK-N, IWORK, INFO1 ) ! LAPACK CALL
  783. CALL SLACPY( 'A', M, N, Z, LDZ, X, LDX ) ! LAPACK CALL
  784. T_OR_N = 'N'
  785. XSCL1 = WORK(N+1)
  786. XSCL2 = WORK(N+2)
  787. IF ( XSCL1 /= XSCL2 ) THEN
  788. ! This is an exceptional situation. If the
  789. ! data matrices are not scaled and the
  790. ! largest singular value of X overflows.
  791. ! In that case SGEJSV can return the SVD
  792. ! in scaled form. The scaling factor can be used
  793. ! to rescale the data (X and Y).
  794. CALL SLASCL( 'G', 0, 0, XSCL1, XSCL2, M, N, Y, LDY, INFO2 )
  795. END IF
  796. END SELECT
  797. !
  798. IF ( INFO1 > 0 ) THEN
  799. ! The SVD selected subroutine did not converge.
  800. ! Return with an error code.
  801. INFO = 2
  802. RETURN
  803. END IF
  804. !
  805. IF ( WORK(1) == ZERO ) THEN
  806. ! The largest computed singular value of (scaled)
  807. ! X is zero. Return error code -8
  808. ! (the 8th input variable had an illegal value).
  809. K = 0
  810. INFO = -8
  811. CALL XERBLA('SGEDMD',-INFO)
  812. RETURN
  813. END IF
  814. !
  815. !<3> Determine the numerical rank of the data
  816. ! snapshots matrix X. This depends on the
  817. ! parameters NRNK and TOL.
  818. SELECT CASE ( NRNK )
  819. CASE ( -1 )
  820. K = 1
  821. DO i = 2, NUMRNK
  822. IF ( ( WORK(i) <= WORK(1)*TOL ) .OR. &
  823. ( WORK(i) <= SMALL ) ) EXIT
  824. K = K + 1
  825. END DO
  826. CASE ( -2 )
  827. K = 1
  828. DO i = 1, NUMRNK-1
  829. IF ( ( WORK(i+1) <= WORK(i)*TOL ) .OR. &
  830. ( WORK(i) <= SMALL ) ) EXIT
  831. K = K + 1
  832. END DO
  833. CASE DEFAULT
  834. K = 1
  835. DO i = 2, NRNK
  836. IF ( WORK(i) <= SMALL ) EXIT
  837. K = K + 1
  838. END DO
  839. END SELECT
  840. ! Now, U = X(1:M,1:K) is the SVD/POD basis for the
  841. ! snapshot data in the input matrix X.
  842. !<4> Compute the Rayleigh quotient S = U^T * A * U.
  843. ! Depending on the requested outputs, the computation
  844. ! is organized to compute additional auxiliary
  845. ! matrices (for the residuals and refinements).
  846. !
  847. ! In all formulas below, we need V_k*Sigma_k^(-1)
  848. ! where either V_k is in W(1:N,1:K), or V_k^T is in
  849. ! W(1:K,1:N). Here Sigma_k=diag(WORK(1:K)).
  850. IF ( LSAME(T_OR_N, 'N') ) THEN
  851. DO i = 1, K
  852. CALL SSCAL( N, ONE/WORK(i), W(1,i), 1 ) ! BLAS CALL
  853. ! W(1:N,i) = (ONE/WORK(i)) * W(1:N,i) ! INTRINSIC
  854. END DO
  855. ELSE
  856. ! This non-unit stride access is due to the fact
  857. ! that SGESVD, SGESVDQ and SGESDD return the
  858. ! transposed matrix of the right singular vectors.
  859. !DO i = 1, K
  860. ! CALL SSCAL( N, ONE/WORK(i), W(i,1), LDW ) ! BLAS CALL
  861. ! ! W(i,1:N) = (ONE/WORK(i)) * W(i,1:N) ! INTRINSIC
  862. !END DO
  863. DO i = 1, K
  864. WORK(N+i) = ONE/WORK(i)
  865. END DO
  866. DO j = 1, N
  867. DO i = 1, K
  868. W(i,j) = (WORK(N+i))*W(i,j)
  869. END DO
  870. END DO
  871. END IF
  872. !
  873. IF ( WNTREF ) THEN
  874. !
  875. ! Need A*U(:,1:K)=Y*V_k*inv(diag(WORK(1:K)))
  876. ! for computing the refined Ritz vectors
  877. ! (optionally, outside SGEDMD).
  878. CALL SGEMM( 'N', T_OR_N, M, K, N, ONE, Y, LDY, W, &
  879. LDW, ZERO, Z, LDZ ) ! BLAS CALL
  880. ! Z(1:M,1:K)=MATMUL(Y(1:M,1:N),TRANSPOSE(W(1:K,1:N))) ! INTRINSIC, for T_OR_N=='T'
  881. ! Z(1:M,1:K)=MATMUL(Y(1:M,1:N),W(1:N,1:K)) ! INTRINSIC, for T_OR_N=='N'
  882. !
  883. ! At this point Z contains
  884. ! A * U(:,1:K) = Y * V_k * Sigma_k^(-1), and
  885. ! this is needed for computing the residuals.
  886. ! This matrix is returned in the array B and
  887. ! it can be used to compute refined Ritz vectors.
  888. CALL SLACPY( 'A', M, K, Z, LDZ, B, LDB ) ! BLAS CALL
  889. ! B(1:M,1:K) = Z(1:M,1:K) ! INTRINSIC
  890. CALL SGEMM( 'T', 'N', K, K, M, ONE, X, LDX, Z, &
  891. LDZ, ZERO, S, LDS ) ! BLAS CALL
  892. ! S(1:K,1:K) = MATMUL(TANSPOSE(X(1:M,1:K)),Z(1:M,1:K)) ! INTRINSIC
  893. ! At this point S = U^T * A * U is the Rayleigh quotient.
  894. ELSE
  895. ! A * U(:,1:K) is not explicitly needed and the
  896. ! computation is organized differently. The Rayleigh
  897. ! quotient is computed more efficiently.
  898. CALL SGEMM( 'T', 'N', K, N, M, ONE, X, LDX, Y, LDY, &
  899. ZERO, Z, LDZ ) ! BLAS CALL
  900. ! Z(1:K,1:N) = MATMUL( TRANSPOSE(X(1:M,1:K)), Y(1:M,1:N) ) ! INTRINSIC
  901. ! In the two SGEMM calls here, can use K for LDZ
  902. CALL SGEMM( 'N', T_OR_N, K, K, N, ONE, Z, LDZ, W, &
  903. LDW, ZERO, S, LDS ) ! BLAS CALL
  904. ! S(1:K,1:K) = MATMUL(Z(1:K,1:N),TRANSPOSE(W(1:K,1:N))) ! INTRINSIC, for T_OR_N=='T'
  905. ! S(1:K,1:K) = MATMUL(Z(1:K,1:N),(W(1:N,1:K))) ! INTRINSIC, for T_OR_N=='N'
  906. ! At this point S = U^T * A * U is the Rayleigh quotient.
  907. ! If the residuals are requested, save scaled V_k into Z.
  908. ! Recall that V_k or V_k^T is stored in W.
  909. IF ( WNTRES .OR. WNTEX ) THEN
  910. IF ( LSAME(T_OR_N, 'N') ) THEN
  911. CALL SLACPY( 'A', N, K, W, LDW, Z, LDZ )
  912. ELSE
  913. CALL SLACPY( 'A', K, N, W, LDW, Z, LDZ )
  914. END IF
  915. END IF
  916. END IF
  917. !
  918. !<5> Compute the Ritz values and (if requested) the
  919. ! right eigenvectors of the Rayleigh quotient.
  920. !
  921. CALL SGEEV( 'N', JOBZL, K, S, LDS, REIG, IMEIG, W, &
  922. LDW, W, LDW, WORK(N+1), LWORK-N, INFO1 ) ! LAPACK CALL
  923. !
  924. ! W(1:K,1:K) contains the eigenvectors of the Rayleigh
  925. ! quotient. Even in the case of complex spectrum, all
  926. ! computation is done in real arithmetic. REIG and
  927. ! IMEIG are the real and the imaginary parts of the
  928. ! eigenvalues, so that the spectrum is given as
  929. ! REIG(:) + sqrt(-1)*IMEIG(:). Complex conjugate pairs
  930. ! are listed at consecutive positions. For such a
  931. ! complex conjugate pair of the eigenvalues, the
  932. ! corresponding eigenvectors are also a complex
  933. ! conjugate pair with the real and imaginary parts
  934. ! stored column-wise in W at the corresponding
  935. ! consecutive column indices. See the description of Z.
  936. ! Also, see the description of SGEEV.
  937. IF ( INFO1 > 0 ) THEN
  938. ! SGEEV failed to compute the eigenvalues and
  939. ! eigenvectors of the Rayleigh quotient.
  940. INFO = 3
  941. RETURN
  942. END IF
  943. !
  944. ! <6> Compute the eigenvectors (if requested) and,
  945. ! the residuals (if requested).
  946. !
  947. IF ( WNTVEC .OR. WNTEX ) THEN
  948. IF ( WNTRES ) THEN
  949. IF ( WNTREF ) THEN
  950. ! Here, if the refinement is requested, we have
  951. ! A*U(:,1:K) already computed and stored in Z.
  952. ! For the residuals, need Y = A * U(:,1;K) * W.
  953. CALL SGEMM( 'N', 'N', M, K, K, ONE, Z, LDZ, W, &
  954. LDW, ZERO, Y, LDY ) ! BLAS CALL
  955. ! Y(1:M,1:K) = Z(1:M,1:K) * W(1:K,1:K) ! INTRINSIC
  956. ! This frees Z; Y contains A * U(:,1:K) * W.
  957. ELSE
  958. ! Compute S = V_k * Sigma_k^(-1) * W, where
  959. ! V_k * Sigma_k^(-1) is stored in Z
  960. CALL SGEMM( T_OR_N, 'N', N, K, K, ONE, Z, LDZ, &
  961. W, LDW, ZERO, S, LDS )
  962. ! Then, compute Z = Y * S =
  963. ! = Y * V_k * Sigma_k^(-1) * W(1:K,1:K) =
  964. ! = A * U(:,1:K) * W(1:K,1:K)
  965. CALL SGEMM( 'N', 'N', M, K, N, ONE, Y, LDY, S, &
  966. LDS, ZERO, Z, LDZ )
  967. ! Save a copy of Z into Y and free Z for holding
  968. ! the Ritz vectors.
  969. CALL SLACPY( 'A', M, K, Z, LDZ, Y, LDY )
  970. IF ( WNTEX ) CALL SLACPY( 'A', M, K, Z, LDZ, B, LDB )
  971. END IF
  972. ELSE IF ( WNTEX ) THEN
  973. ! Compute S = V_k * Sigma_k^(-1) * W, where
  974. ! V_k * Sigma_k^(-1) is stored in Z
  975. CALL SGEMM( T_OR_N, 'N', N, K, K, ONE, Z, LDZ, &
  976. W, LDW, ZERO, S, LDS )
  977. ! Then, compute Z = Y * S =
  978. ! = Y * V_k * Sigma_k^(-1) * W(1:K,1:K) =
  979. ! = A * U(:,1:K) * W(1:K,1:K)
  980. CALL SGEMM( 'N', 'N', M, K, N, ONE, Y, LDY, S, &
  981. LDS, ZERO, B, LDB )
  982. ! The above call replaces the following two calls
  983. ! that were used in the developing-testing phase.
  984. ! CALL SGEMM( 'N', 'N', M, K, N, ONE, Y, LDY, S, &
  985. ! LDS, ZERO, Z, LDZ)
  986. ! Save a copy of Z into B and free Z for holding
  987. ! the Ritz vectors.
  988. ! CALL SLACPY( 'A', M, K, Z, LDZ, B, LDB )
  989. END IF
  990. !
  991. ! Compute the real form of the Ritz vectors
  992. IF ( WNTVEC ) CALL SGEMM( 'N', 'N', M, K, K, ONE, X, LDX, W, LDW, &
  993. ZERO, Z, LDZ ) ! BLAS CALL
  994. ! Z(1:M,1:K) = MATMUL(X(1:M,1:K), W(1:K,1:K)) ! INTRINSIC
  995. !
  996. IF ( WNTRES ) THEN
  997. i = 1
  998. DO WHILE ( i <= K )
  999. IF ( IMEIG(i) == ZERO ) THEN
  1000. ! have a real eigenvalue with real eigenvector
  1001. CALL SAXPY( M, -REIG(i), Z(1,i), 1, Y(1,i), 1 ) ! BLAS CALL
  1002. ! Y(1:M,i) = Y(1:M,i) - REIG(i) * Z(1:M,i) ! INTRINSIC
  1003. RES(i) = SNRM2( M, Y(1,i), 1 ) ! BLAS CALL
  1004. i = i + 1
  1005. ELSE
  1006. ! Have a complex conjugate pair
  1007. ! REIG(i) +- sqrt(-1)*IMEIG(i).
  1008. ! Since all computation is done in real
  1009. ! arithmetic, the formula for the residual
  1010. ! is recast for real representation of the
  1011. ! complex conjugate eigenpair. See the
  1012. ! description of RES.
  1013. AB(1,1) = REIG(i)
  1014. AB(2,1) = -IMEIG(i)
  1015. AB(1,2) = IMEIG(i)
  1016. AB(2,2) = REIG(i)
  1017. CALL SGEMM( 'N', 'N', M, 2, 2, -ONE, Z(1,i), &
  1018. LDZ, AB, 2, ONE, Y(1,i), LDY ) ! BLAS CALL
  1019. ! Y(1:M,i:i+1) = Y(1:M,i:i+1) - Z(1:M,i:i+1) * AB ! INTRINSIC
  1020. RES(i) = SLANGE( 'F', M, 2, Y(1,i), LDY, &
  1021. WORK(N+1) ) ! LAPACK CALL
  1022. RES(i+1) = RES(i)
  1023. i = i + 2
  1024. END IF
  1025. END DO
  1026. END IF
  1027. END IF
  1028. !
  1029. IF ( WHTSVD == 4 ) THEN
  1030. WORK(N+1) = XSCL1
  1031. WORK(N+2) = XSCL2
  1032. END IF
  1033. !
  1034. ! Successful exit.
  1035. IF ( .NOT. BADXY ) THEN
  1036. INFO = 0
  1037. ELSE
  1038. ! A warning on possible data inconsistency.
  1039. ! This should be a rare event.
  1040. INFO = 4
  1041. END IF
  1042. !............................................................
  1043. RETURN
  1044. ! ......
  1045. END SUBROUTINE SGEDMD