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.

schkdmd.f90 29 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792
  1. ! This is a test program for checking the implementations of
  2. ! the implementations of the following subroutines
  3. !
  4. ! SGEDMD for computation of the
  5. ! Dynamic Mode Decomposition (DMD)
  6. ! SGEDMDQ for computation of a
  7. ! QR factorization based compressed DMD
  8. !
  9. ! Developed and supported by:
  10. ! ===========================
  11. ! Developed and coded by Zlatko Drmac, Faculty of Science,
  12. ! University of Zagreb; drmac@math.hr
  13. ! In cooperation with
  14. ! AIMdyn Inc., Santa Barbara, CA.
  15. ! ========================================================
  16. ! How to run the code (compiler, link info)
  17. ! ========================================================
  18. ! Compile as FORTRAN 90 (or later) and link with BLAS and
  19. ! LAPACK libraries.
  20. ! NOTE: The code is developed and tested on top of the
  21. ! Intel MKL library (versions 2022.0.3 and 2022.2.0),
  22. ! using the Intel Fortran compiler.
  23. !
  24. ! For developers of the C++ implementation
  25. ! ========================================================
  26. ! See the LAPACK++ and Template Numerical Toolkit (TNT)
  27. !
  28. ! Note on a development of the GPU HP implementation
  29. ! ========================================================
  30. ! Work in progress. See CUDA, MAGMA, SLATE.
  31. ! NOTE: The four SVD subroutines used in this code are
  32. ! included as a part of R&D and for the completeness.
  33. ! This was also an opportunity to test those SVD codes.
  34. ! If the scaling option is used all four are essentially
  35. ! equally good. For implementations on HP platforms,
  36. ! one can use whichever SVD is available.
  37. !... .........................................................
  38. ! NOTE:
  39. ! When using the Intel MKL 2022.0.3 the subroutine xGESVDQ
  40. ! (optionally used in xGEDMD) may cause access violation
  41. ! error for x = S, D, C, Z, but only if called with the
  42. ! work space query. (At least in our Windows 10 MSVS 2019.)
  43. ! The problem can be mitigated by downloading the source
  44. ! code of xGESVDQ from the LAPACK repository and use it
  45. ! localy instead of the one in the MKL. This seems to
  46. ! indicate that the problem is indeed in the MKL.
  47. ! This problem did not appear whith Intel MKL 2022.2.0.
  48. !
  49. ! NOTE:
  50. ! xGESDD seems to have a problem with workspace. In some
  51. ! cases the length of the optimal workspace is returned
  52. ! smaller than the minimal workspace, as specified in the
  53. ! code. As a precaution, all optimal workspaces are
  54. ! set as MAX(minimal, optimal).
  55. ! Latest implementations of complex xGESDD have different
  56. ! length of the real worksapce. We use max value over
  57. ! two versions.
  58. !............................................................
  59. !............................................................
  60. !
  61. PROGRAM DMD_TEST
  62. use iso_fortran_env, only: real32
  63. IMPLICIT NONE
  64. integer, parameter :: WP = real32
  65. !............................................................
  66. REAL(KIND=WP), PARAMETER :: ONE = 1.0_WP
  67. REAL(KIND=WP), PARAMETER :: ZERO = 0.0_WP
  68. !............................................................
  69. REAL(KIND=WP), ALLOCATABLE, DIMENSION(:,:) :: &
  70. A, AC, EIGA, LAMBDA, LAMBDAQ, F, F1, F2,&
  71. Z, Z1, S, AU, W, VA, X, X0, Y, Y0, Y1
  72. REAL(KIND=WP), ALLOCATABLE, DIMENSION(:) :: &
  73. DA, DL, DR, REIG, REIGA, REIGQ, IEIG, &
  74. IEIGA, IEIGQ, RES, RES1, RESEX, SINGVX,&
  75. SINGVQX, WORK
  76. INTEGER , ALLOCATABLE, DIMENSION(:) :: IWORK
  77. REAL(KIND=WP) :: AB(2,2), WDUMMY(2)
  78. INTEGER :: IDUMMY(2), ISEED(4), RJOBDATA(8)
  79. REAL(KIND=WP) :: ANORM, COND, CONDL, CONDR, DMAX, EPS, &
  80. TOL, TOL2, SVDIFF, TMP, TMP_AU, &
  81. TMP_FQR, TMP_REZ, TMP_REZQ, TMP_ZXW, &
  82. TMP_EX, XNORM, YNORM
  83. !............................................................
  84. INTEGER :: K, KQ, LDF, LDS, LDA, LDAU, LDW, LDX, LDY, &
  85. LDZ, LIWORK, LWORK, M, N, L, LLOOP, NRNK
  86. INTEGER :: i, iJOBREF, iJOBZ, iSCALE, INFO, KDIFF, &
  87. NFAIL, NFAIL_AU, NFAIL_F_QR, NFAIL_REZ, &
  88. NFAIL_REZQ, NFAIL_SVDIFF, NFAIL_TOTAL, NFAILQ_TOTAL, &
  89. NFAIL_Z_XV, MODE, MODEL, MODER, WHTSVD
  90. INTEGER iNRNK, iWHTSVD, K_TRAJ, LWMINOPT
  91. CHARACTER(LEN=1) GRADE, JOBREF, JOBZ, PIVTNG, RSIGN, &
  92. SCALE, RESIDS, WANTQ, WANTR
  93. LOGICAL TEST_QRDMD
  94. !..... external subroutines (BLAS and LAPACK)
  95. EXTERNAL SAXPY, SGEEV, SGEMM, SGEMV, SLACPY, SLASCL
  96. EXTERNAL SLARNV, SLATMR
  97. !.....external subroutines DMD package, part 1
  98. ! subroutines under test
  99. EXTERNAL SGEDMD, SGEDMDQ
  100. !..... external functions (BLAS and LAPACK)
  101. EXTERNAL SLAMCH, SLANGE, SNRM2
  102. REAL(KIND=WP) :: SLAMCH, SLANGE, SNRM2
  103. EXTERNAL LSAME
  104. LOGICAL LSAME
  105. INTRINSIC ABS, INT, MIN, MAX
  106. !............................................................
  107. ! The test is always in pairs : ( SGEDMD and SGEDMDQ )
  108. ! because the test includes comparing the results (in pairs).
  109. !.....................................................................................
  110. TEST_QRDMD = .TRUE. ! This code by default performs tests on SGEDMDQ
  111. ! Since the QR factorizations based algorithm is designed for
  112. ! single trajectory data, only single trajectory tests will
  113. ! be performed with xGEDMDQ.
  114. WANTQ = 'Q'
  115. WANTR = 'R'
  116. !.................................................................................
  117. EPS = SLAMCH( 'P' ) ! machine precision SP
  118. ! Global counters of failures of some particular tests
  119. NFAIL = 0
  120. NFAIL_REZ = 0
  121. NFAIL_REZQ = 0
  122. NFAIL_Z_XV = 0
  123. NFAIL_F_QR = 0
  124. NFAIL_AU = 0
  125. KDIFF = 0
  126. NFAIL_SVDIFF = 0
  127. NFAIL_TOTAL = 0
  128. NFAILQ_TOTAL = 0
  129. DO LLOOP = 1, 4
  130. WRITE(*,*) 'L Loop Index = ', LLOOP
  131. ! Set the dimensions of the problem ...
  132. WRITE(*,*) 'M = '
  133. READ(*,*) M
  134. WRITE(*,*) M
  135. ! ... and the number of snapshots.
  136. WRITE(*,*) 'N = '
  137. READ(*,*) N
  138. WRITE(*,*) N
  139. ! ... Test the dimensions
  140. IF ( ( MIN(M,N) == 0 ) .OR. ( M < N ) ) THEN
  141. WRITE(*,*) 'Bad dimensions. Required: M >= N > 0.'
  142. STOP
  143. END IF
  144. !.............
  145. ! The seed inside the LLOOP so that each pass can be reproduced easily.
  146. ISEED(1) = 4
  147. ISEED(2) = 3
  148. ISEED(3) = 2
  149. ISEED(4) = 1
  150. LDA = M
  151. LDF = M
  152. LDX = MAX(M,N+1)
  153. LDY = MAX(M,N+1)
  154. LDW = N
  155. LDZ = M
  156. LDAU = MAX(M,N+1)
  157. LDS = N
  158. TMP_ZXW = ZERO
  159. TMP_AU = ZERO
  160. TMP_REZ = ZERO
  161. TMP_REZQ = ZERO
  162. SVDIFF = ZERO
  163. TMP_EX = ZERO
  164. !
  165. ! Test the subroutines on real data snapshots. All
  166. ! computation is done in real arithmetic, even when
  167. ! Koopman eigenvalues and modes are real.
  168. !
  169. ! Allocate memory space
  170. ALLOCATE( A(LDA,M) )
  171. ALLOCATE( AC(LDA,M) )
  172. ALLOCATE( DA(M) )
  173. ALLOCATE( DL(M) )
  174. ALLOCATE( F(LDF,N+1) )
  175. ALLOCATE( F1(LDF,N+1) )
  176. ALLOCATE( F2(LDF,N+1) )
  177. ALLOCATE( X(LDX,N) )
  178. ALLOCATE( X0(LDX,N) )
  179. ALLOCATE( SINGVX(N) )
  180. ALLOCATE( SINGVQX(N) )
  181. ALLOCATE( Y(LDY,N+1) )
  182. ALLOCATE( Y0(LDY,N+1) )
  183. ALLOCATE( Y1(M,N+1) )
  184. ALLOCATE( Z(LDZ,N) )
  185. ALLOCATE( Z1(LDZ,N) )
  186. ALLOCATE( RES(N) )
  187. ALLOCATE( RES1(N) )
  188. ALLOCATE( RESEX(N) )
  189. ALLOCATE( REIG(N) )
  190. ALLOCATE( IEIG(N) )
  191. ALLOCATE( REIGQ(N) )
  192. ALLOCATE( IEIGQ(N) )
  193. ALLOCATE( REIGA(M) )
  194. ALLOCATE( IEIGA(M) )
  195. ALLOCATE( VA(LDA,M) )
  196. ALLOCATE( LAMBDA(N,2) )
  197. ALLOCATE( LAMBDAQ(N,2) )
  198. ALLOCATE( EIGA(M,2) )
  199. ALLOCATE( W(LDW,N) )
  200. ALLOCATE( AU(LDAU,N) )
  201. ALLOCATE( S(N,N) )
  202. TOL = M*EPS
  203. ! This mimics O(M*N)*EPS bound for accumulated roundoff error.
  204. ! The factor 10 is somewhat arbitrary.
  205. TOL2 = 10*M*N*EPS
  206. !.............
  207. DO K_TRAJ = 1, 2
  208. ! Number of intial conditions in the simulation/trajectories (1 or 2)
  209. COND = 1.0D8
  210. DMAX = 1.0D2
  211. RSIGN = 'F'
  212. GRADE = 'N'
  213. MODEL = 6
  214. CONDL = 1.0D2
  215. MODER = 6
  216. CONDR = 1.0D2
  217. PIVTNG = 'N'
  218. ! Loop over all parameter MODE values for ZLATMR (+1,..,+6)
  219. DO MODE = 1, 6
  220. ALLOCATE( IWORK(2*M) )
  221. ALLOCATE(DR(N))
  222. CALL SLATMR( M, M, 'S', ISEED, 'N', DA, MODE, COND, &
  223. DMAX, RSIGN, GRADE, DL, MODEL, CONDL, &
  224. DR, MODER, CONDR, PIVTNG, IWORK, M, M, &
  225. ZERO, -ONE, 'N', A, LDA, IWORK(M+1), INFO )
  226. DEALLOCATE(IWORK)
  227. DEALLOCATE(DR)
  228. LWORK = 4*M+1
  229. ALLOCATE(WORK(LWORK))
  230. AC = A
  231. CALL SGEEV( 'N','V', M, AC, M, REIGA, IEIGA, VA, M, &
  232. VA, M, WORK, LWORK, INFO ) ! LAPACK CALL
  233. DEALLOCATE(WORK)
  234. TMP = ZERO
  235. DO i = 1, M
  236. EIGA(i,1) = REIGA(i)
  237. EIGA(i,2) = IEIGA(i)
  238. TMP = MAX( TMP, SQRT(REIGA(i)**2+IEIGA(i)**2))
  239. END DO
  240. ! Scale A to have the desirable spectral radius.
  241. CALL SLASCL( 'G', 0, 0, TMP, ONE, M, M, A, M, INFO )
  242. CALL SLASCL( 'G', 0, 0, TMP, ONE, M, 2, EIGA, M, INFO )
  243. ! Compute the norm of A
  244. ANORM = SLANGE( 'F', N, N, A, M, WDUMMY )
  245. IF ( K_TRAJ == 2 ) THEN
  246. ! generate data with two inital conditions
  247. CALL SLARNV(2, ISEED, M, F1(1,1) )
  248. F1(1:M,1) = 1.0E-10*F1(1:M,1)
  249. DO i = 1, N/2
  250. CALL SGEMV( 'N', M, M, ONE, A, M, F1(1,i), 1, ZERO, &
  251. F1(1,i+1), 1 )
  252. END DO
  253. X0(1:M,1:N/2) = F1(1:M,1:N/2)
  254. Y0(1:M,1:N/2) = F1(1:M,2:N/2+1)
  255. CALL SLARNV(2, ISEED, M, F1(1,1) )
  256. DO i = 1, N-N/2
  257. CALL SGEMV( 'N', M, M, ONE, A, M, F1(1,i), 1, ZERO, &
  258. F1(1,i+1), 1 )
  259. END DO
  260. X0(1:M,N/2+1:N) = F1(1:M,1:N-N/2)
  261. Y0(1:M,N/2+1:N) = F1(1:M,2:N-N/2+1)
  262. ELSE
  263. ! single trajectory
  264. CALL SLARNV(2, ISEED, M, F(1,1) )
  265. DO i = 1, N
  266. CALL SGEMV( 'N', M, M, ONE, A, M, F(1,i), 1, ZERO, &
  267. F(1,i+1), 1 )
  268. END DO
  269. X0(1:M,1:N) = F(1:M,1:N)
  270. Y0(1:M,1:N) = F(1:M,2:N+1)
  271. END IF
  272. XNORM = SLANGE( 'F', M, N, X0, LDX, WDUMMY )
  273. YNORM = SLANGE( 'F', M, N, Y0, LDX, WDUMMY )
  274. !............................................................
  275. DO iJOBZ = 1, 4
  276. SELECT CASE ( iJOBZ )
  277. CASE(1)
  278. JOBZ = 'V' ! Ritz vectors will be computed
  279. RESIDS = 'R' ! Residuals will be computed
  280. CASE(2)
  281. JOBZ = 'V'
  282. RESIDS = 'N'
  283. CASE(3)
  284. JOBZ = 'F' ! Ritz vectors in factored form
  285. RESIDS = 'N'
  286. CASE(4)
  287. JOBZ = 'N'
  288. RESIDS = 'N'
  289. END SELECT
  290. DO iJOBREF = 1, 3
  291. SELECT CASE ( iJOBREF )
  292. CASE(1)
  293. JOBREF = 'R' ! Data for refined Ritz vectors
  294. CASE(2)
  295. JOBREF = 'E' ! Exact DMD vectors
  296. CASE(3)
  297. JOBREF = 'N'
  298. END SELECT
  299. DO iSCALE = 1, 4
  300. SELECT CASE ( iSCALE )
  301. CASE(1)
  302. SCALE = 'S' ! X data normalized
  303. CASE(2)
  304. SCALE = 'C' ! X normalized, consist. check
  305. CASE(3)
  306. SCALE = 'Y' ! Y data normalized
  307. CASE(4)
  308. SCALE = 'N'
  309. END SELECT
  310. DO iNRNK = -1, -2, -1
  311. ! Two truncation strategies. The "-2" case for R&D
  312. ! purposes only - it uses possibly low accuracy small
  313. ! singular values, in which case the formulas used in
  314. ! the DMD are highly sensitive.
  315. NRNK = iNRNK
  316. DO iWHTSVD = 1, 4
  317. ! Check all four options to compute the POD basis
  318. ! via the SVD.
  319. WHTSVD = iWHTSVD
  320. DO LWMINOPT = 1, 2
  321. ! Workspace query for the minimal (1) and for the optimal
  322. ! (2) workspace lengths determined by workspace query.
  323. X(1:M,1:N) = X0(1:M,1:N)
  324. Y(1:M,1:N) = Y0(1:M,1:N)
  325. ! SGEDMD: Workspace query and workspace allocation
  326. CALL SGEDMD( SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, M, &
  327. N, X, LDX, Y, LDY, NRNK, TOL, K, REIG, IEIG, Z, &
  328. LDZ, RES, AU, LDAU, W, LDW, S, LDS, WDUMMY, -1, &
  329. IDUMMY, -1, INFO )
  330. LIWORK = IDUMMY(1)
  331. ALLOCATE( IWORK(LIWORK) )
  332. LWORK = INT(WDUMMY(LWMINOPT))
  333. ALLOCATE( WORK(LWORK) )
  334. ! SGEDMD test: CALL SGEDMD
  335. CALL SGEDMD( SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, M, &
  336. N, X, LDX, Y, LDY, NRNK, TOL, K, REIG, IEIG, Z, &
  337. LDZ, RES, AU, LDAU, W, LDW, S, LDS, WORK, LWORK,&
  338. IWORK, LIWORK, INFO )
  339. SINGVX(1:N) = WORK(1:N)
  340. !...... SGEDMD check point
  341. IF ( LSAME(JOBZ,'V') ) THEN
  342. ! Check that Z = X*W, on return from SGEDMD
  343. ! This checks that the returned aigenvectors in Z are
  344. ! the product of the SVD'POD basis returned in X
  345. ! and the eigenvectors of the rayleigh quotient
  346. ! returned in W
  347. CALL SGEMM( 'N', 'N', M, K, K, ONE, X, LDX, W, LDW, &
  348. ZERO, Z1, LDZ )
  349. TMP = ZERO
  350. DO i = 1, K
  351. CALL SAXPY( M, -ONE, Z(1,i), 1, Z1(1,i), 1)
  352. TMP = MAX(TMP, SNRM2( M, Z1(1,i), 1 ) )
  353. END DO
  354. TMP_ZXW = MAX(TMP_ZXW, TMP )
  355. IF ( TMP_ZXW > 10*M*EPS ) THEN
  356. NFAIL_Z_XV = NFAIL_Z_XV + 1
  357. END IF
  358. END IF
  359. !...... SGEDMD check point
  360. IF ( LSAME(JOBREF,'R') ) THEN
  361. ! The matrix A*U is returned for computing refined Ritz vectors.
  362. ! Check that A*U is computed correctly using the formula
  363. ! A*U = Y * V * inv(SIGMA). This depends on the
  364. ! accuracy in the computed singular values and vectors of X.
  365. ! See the paper for an error analysis.
  366. ! Note that the left singular vectors of the input matrix X
  367. ! are returned in the array X.
  368. CALL SGEMM( 'N', 'N', M, K, M, ONE, A, LDA, X, LDX, &
  369. ZERO, Z1, LDZ )
  370. TMP = ZERO
  371. DO i = 1, K
  372. CALL SAXPY( M, -ONE, AU(1,i), 1, Z1(1,i), 1)
  373. TMP = MAX( TMP, SNRM2( M, Z1(1,i),1 ) * &
  374. SINGVX(K)/(ANORM*SINGVX(1)) )
  375. END DO
  376. TMP_AU = MAX( TMP_AU, TMP )
  377. IF ( TMP > TOL2 ) THEN
  378. NFAIL_AU = NFAIL_AU + 1
  379. END IF
  380. ELSEIF ( LSAME(JOBREF,'E') ) THEN
  381. ! The unscaled vectors of the Exact DMD are computed.
  382. ! This option is included for the sake of completeness,
  383. ! for users who prefer the Exact DMD vectors. The
  384. ! returned vectors are in the real form, in the same way
  385. ! as the Ritz vectors. Here we just save the vectors
  386. ! and test them separately using a Matlab script.
  387. CALL SGEMM( 'N', 'N', M, K, M, ONE, A, LDA, AU, LDAU, ZERO, Y1, M )
  388. i=1
  389. DO WHILE ( i <= K )
  390. IF ( IEIG(i) == ZERO ) THEN
  391. ! have a real eigenvalue with real eigenvector
  392. CALL SAXPY( M, -REIG(i), AU(1,i), 1, Y1(1,i), 1 )
  393. RESEX(i) = SNRM2( M, Y1(1,i), 1) / SNRM2(M,AU(1,i),1)
  394. i = i + 1
  395. ELSE
  396. ! Have a complex conjugate pair
  397. ! REIG(i) +- sqrt(-1)*IMEIG(i).
  398. ! Since all computation is done in real
  399. ! arithmetic, the formula for the residual
  400. ! is recast for real representation of the
  401. ! complex conjugate eigenpair. See the
  402. ! description of RES.
  403. AB(1,1) = REIG(i)
  404. AB(2,1) = -IEIG(i)
  405. AB(1,2) = IEIG(i)
  406. AB(2,2) = REIG(i)
  407. CALL SGEMM( 'N', 'N', M, 2, 2, -ONE, AU(1,i), &
  408. M, AB, 2, ONE, Y1(1,i), M )
  409. RESEX(i) = SLANGE( 'F', M, 2, Y1(1,i), M, &
  410. WORK )/ SLANGE( 'F', M, 2, AU(1,i), M, &
  411. WORK )
  412. RESEX(i+1) = RESEX(i)
  413. i = i + 2
  414. END IF
  415. END DO
  416. END IF
  417. !...... SGEDMD check point
  418. IF ( LSAME(RESIDS, 'R') ) THEN
  419. ! Compare the residuals returned by SGEDMD with the
  420. ! explicitly computed residuals using the matrix A.
  421. ! Compute explicitly Y1 = A*Z
  422. CALL SGEMM( 'N', 'N', M, K, M, ONE, A, LDA, Z, LDZ, ZERO, Y1, M )
  423. ! ... and then A*Z(:,i) - LAMBDA(i)*Z(:,i), using the real forms
  424. ! of the invariant subspaces that correspond to complex conjugate
  425. ! pairs of eigencalues. (See the description of Z in SGEDMD,)
  426. i = 1
  427. DO WHILE ( i <= K )
  428. IF ( IEIG(i) == ZERO ) THEN
  429. ! have a real eigenvalue with real eigenvector
  430. CALL SAXPY( M, -REIG(i), Z(1,i), 1, Y1(1,i), 1 )
  431. RES1(i) = SNRM2( M, Y1(1,i), 1)
  432. i = i + 1
  433. ELSE
  434. ! Have a complex conjugate pair
  435. ! REIG(i) +- sqrt(-1)*IMEIG(i).
  436. ! Since all computation is done in real
  437. ! arithmetic, the formula for the residual
  438. ! is recast for real representation of the
  439. ! complex conjugate eigenpair. See the
  440. ! description of RES.
  441. AB(1,1) = REIG(i)
  442. AB(2,1) = -IEIG(i)
  443. AB(1,2) = IEIG(i)
  444. AB(2,2) = REIG(i)
  445. CALL SGEMM( 'N', 'N', M, 2, 2, -ONE, Z(1,i), &
  446. M, AB, 2, ONE, Y1(1,i), M )
  447. RES1(i) = SLANGE( 'F', M, 2, Y1(1,i), M, &
  448. WORK )
  449. RES1(i+1) = RES1(i)
  450. i = i + 2
  451. END IF
  452. END DO
  453. TMP = ZERO
  454. DO i = 1, K
  455. TMP = MAX( TMP, ABS(RES(i) - RES1(i)) * &
  456. SINGVX(K)/(ANORM*SINGVX(1)) )
  457. END DO
  458. TMP_REZ = MAX( TMP_REZ, TMP )
  459. IF ( TMP > TOL2 ) THEN
  460. NFAIL_REZ = NFAIL_REZ + 1
  461. END IF
  462. IF ( LSAME(JOBREF,'E') ) THEN
  463. TMP = ZERO
  464. DO i = 1, K
  465. TMP = MAX( TMP, ABS(RES1(i) - RESEX(i))/(RES1(i)+RESEX(i)) )
  466. END DO
  467. TMP_EX = MAX(TMP_EX,TMP)
  468. END IF
  469. END IF
  470. ! ... store the results for inspection
  471. DO i = 1, K
  472. LAMBDA(i,1) = REIG(i)
  473. LAMBDA(i,2) = IEIG(i)
  474. END DO
  475. DEALLOCATE(IWORK)
  476. DEALLOCATE(WORK)
  477. !======================================================================
  478. ! Now test the SGEDMDQ, if requested.
  479. !======================================================================
  480. IF ( TEST_QRDMD .AND. (K_TRAJ == 1) ) THEN
  481. RJOBDATA(2) = 1
  482. F1 = F
  483. ! SGEDMDQ test: Workspace query and workspace allocation
  484. CALL SGEDMDQ( SCALE, JOBZ, RESIDS, WANTQ, WANTR, &
  485. JOBREF, WHTSVD, M, N+1, F1, LDF, X, LDX, Y, &
  486. LDY, NRNK, TOL, KQ, REIGQ, IEIGQ, Z, LDZ, &
  487. RES, AU, LDAU, W, LDW, S, LDS, WDUMMY, &
  488. -1, IDUMMY, -1, INFO )
  489. LIWORK = IDUMMY(1)
  490. ALLOCATE( IWORK(LIWORK) )
  491. LWORK = INT(WDUMMY(LWMINOPT))
  492. ALLOCATE(WORK(LWORK))
  493. ! SGEDMDQ test: CALL SGEDMDQ
  494. CALL SGEDMDQ( SCALE, JOBZ, RESIDS, WANTQ, WANTR, &
  495. JOBREF, WHTSVD, M, N+1, F1, LDF, X, LDX, Y, &
  496. LDY, NRNK, TOL, KQ, REIGQ, IEIGQ, Z, LDZ, &
  497. RES, AU, LDAU, W, LDW, S, LDS, &
  498. WORK, LWORK, IWORK, LIWORK, INFO )
  499. SINGVQX(1:KQ) = WORK(MIN(M,N+1)+1: MIN(M,N+1)+KQ)
  500. !..... SGEDMDQ check point
  501. IF ( KQ /= K ) THEN
  502. KDIFF = KDIFF+1
  503. END IF
  504. TMP = ZERO
  505. DO i = 1, MIN(K, KQ)
  506. TMP = MAX(TMP, ABS(SINGVX(i)-SINGVQX(i)) / &
  507. SINGVX(1) )
  508. END DO
  509. SVDIFF = MAX( SVDIFF, TMP )
  510. IF ( TMP > M*N*EPS ) THEN
  511. NFAIL_SVDIFF = NFAIL_SVDIFF + 1
  512. END IF
  513. !..... SGEDMDQ check point
  514. IF ( LSAME(WANTQ,'Q') .AND. LSAME(WANTR,'R') ) THEN
  515. ! Check that the QR factors are computed and returned
  516. ! as requested. The residual ||F-Q*R||_F / ||F||_F
  517. ! is compared to M*N*EPS.
  518. F2 = F
  519. CALL SGEMM( 'N', 'N', M, N+1, MIN(M,N+1), -ONE, F1, &
  520. LDF, Y, LDY, ONE, F2, LDF )
  521. TMP_FQR = SLANGE( 'F', M, N+1, F2, LDF, WORK ) / &
  522. SLANGE( 'F', M, N+1, F, LDF, WORK )
  523. IF ( TMP_FQR > TOL2 ) THEN
  524. NFAIL_F_QR = NFAIL_F_QR + 1
  525. END IF
  526. END IF
  527. !..... SGEDMDQ checkpoint
  528. IF ( LSAME(RESIDS, 'R') ) THEN
  529. ! Compare the residuals returned by SGEDMDQ with the
  530. ! explicitly computed residuals using the matrix A.
  531. ! Compute explicitly Y1 = A*Z
  532. CALL SGEMM( 'N', 'N', M, KQ, M, ONE, A, M, Z, M, ZERO, Y1, M )
  533. ! ... and then A*Z(:,i) - LAMBDA(i)*Z(:,i), using the real forms
  534. ! of the invariant subspaces that correspond to complex conjugate
  535. ! pairs of eigencalues. (See the description of Z in SGEDMDQ)
  536. i = 1
  537. DO WHILE ( i <= KQ )
  538. IF ( IEIGQ(i) == ZERO ) THEN
  539. ! have a real eigenvalue with real eigenvector
  540. CALL SAXPY( M, -REIGQ(i), Z(1,i), 1, Y1(1,i), 1 )
  541. ! Y(1:M,i) = Y(1:M,i) - REIG(i)*Z(1:M,i)
  542. RES1(i) = SNRM2( M, Y1(1,i), 1)
  543. i = i + 1
  544. ELSE
  545. ! Have a complex conjugate pair
  546. ! REIG(i) +- sqrt(-1)*IMEIG(i).
  547. ! Since all computation is done in real
  548. ! arithmetic, the formula for the residual
  549. ! is recast for real representation of the
  550. ! complex conjugate eigenpair. See the
  551. ! description of RES.
  552. AB(1,1) = REIGQ(i)
  553. AB(2,1) = -IEIGQ(i)
  554. AB(1,2) = IEIGQ(i)
  555. AB(2,2) = REIGQ(i)
  556. CALL SGEMM( 'N', 'N', M, 2, 2, -ONE, Z(1,i), &
  557. M, AB, 2, ONE, Y1(1,i), M ) ! BLAS CALL
  558. ! Y(1:M,i:i+1) = Y(1:M,i:i+1) - Z(1:M,i:i+1) * AB ! INTRINSIC
  559. RES1(i) = SLANGE( 'F', M, 2, Y1(1,i), M, &
  560. WORK ) ! LAPACK CALL
  561. RES1(i+1) = RES1(i)
  562. i = i + 2
  563. END IF
  564. END DO
  565. TMP = ZERO
  566. DO i = 1, KQ
  567. TMP = MAX( TMP, ABS(RES(i) - RES1(i)) * &
  568. SINGVQX(K)/(ANORM*SINGVQX(1)) )
  569. END DO
  570. TMP_REZQ = MAX( TMP_REZQ, TMP )
  571. IF ( TMP > TOL2 ) THEN
  572. NFAIL_REZQ = NFAIL_REZQ + 1
  573. END IF
  574. END IF
  575. DO i = 1, KQ
  576. LAMBDAQ(i,1) = REIGQ(i)
  577. LAMBDAQ(i,2) = IEIGQ(i)
  578. END DO
  579. DEALLOCATE(WORK)
  580. DEALLOCATE(IWORK)
  581. END IF ! TEST_QRDMD
  582. !======================================================================
  583. END DO ! LWMINOPT
  584. !write(*,*) 'LWMINOPT loop completed'
  585. END DO ! WHTSVD LOOP
  586. !write(*,*) 'WHTSVD loop completed'
  587. END DO ! NRNK LOOP
  588. !write(*,*) 'NRNK loop completed'
  589. END DO ! SCALE LOOP
  590. !write(*,*) 'SCALE loop completed'
  591. END DO ! JOBF LOOP
  592. !write(*,*) 'JOBREF loop completed'
  593. END DO ! JOBZ LOOP
  594. !write(*,*) 'JOBZ loop completed'
  595. END DO ! MODE -6:6
  596. !write(*,*) 'MODE loop completed'
  597. END DO ! 1 or 2 trajectories
  598. !write(*,*) 'trajectories loop completed'
  599. DEALLOCATE(A)
  600. DEALLOCATE(AC)
  601. DEALLOCATE(DA)
  602. DEALLOCATE(DL)
  603. DEALLOCATE(F)
  604. DEALLOCATE(F1)
  605. DEALLOCATE(F2)
  606. DEALLOCATE(X)
  607. DEALLOCATE(X0)
  608. DEALLOCATE(SINGVX)
  609. DEALLOCATE(SINGVQX)
  610. DEALLOCATE(Y)
  611. DEALLOCATE(Y0)
  612. DEALLOCATE(Y1)
  613. DEALLOCATE(Z)
  614. DEALLOCATE(Z1)
  615. DEALLOCATE(RES)
  616. DEALLOCATE(RES1)
  617. DEALLOCATE(RESEX)
  618. DEALLOCATE(REIG)
  619. DEALLOCATE(IEIG)
  620. DEALLOCATE(REIGQ)
  621. DEALLOCATE(IEIGQ)
  622. DEALLOCATE(REIGA)
  623. DEALLOCATE(IEIGA)
  624. DEALLOCATE(VA)
  625. DEALLOCATE(LAMBDA)
  626. DEALLOCATE(LAMBDAQ)
  627. DEALLOCATE(EIGA)
  628. DEALLOCATE(W)
  629. DEALLOCATE(AU)
  630. DEALLOCATE(S)
  631. !............................................................
  632. ! Generate random M-by-M matrix A. Use DLATMR from
  633. END DO ! LLOOP
  634. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  635. WRITE(*,*) ' Test summary for SGEDMD :'
  636. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  637. WRITE(*,*)
  638. IF ( NFAIL_Z_XV == 0 ) THEN
  639. WRITE(*,*) '>>>> Z - U*V test PASSED.'
  640. ELSE
  641. WRITE(*,*) 'Z - U*V test FAILED ', NFAIL_Z_XV, ' time(s)'
  642. WRITE(*,*) 'Max error ||Z-U*V||_F was ', TMP_ZXW
  643. NFAIL_TOTAL = NFAIL_TOTAL + NFAIL_Z_XV
  644. END IF
  645. IF ( NFAIL_AU == 0 ) THEN
  646. WRITE(*,*) '>>>> A*U test PASSED. '
  647. ELSE
  648. WRITE(*,*) 'A*U test FAILED ', NFAIL_AU, ' time(s)'
  649. WRITE(*,*) 'Max A*U test adjusted error measure was ', TMP_AU
  650. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  651. NFAIL_TOTAL = NFAIL_TOTAL + NFAIL_AU
  652. END IF
  653. IF ( NFAIL_REZ == 0 ) THEN
  654. WRITE(*,*) '>>>> Rezidual computation test PASSED.'
  655. ELSE
  656. WRITE(*,*) 'Rezidual computation test FAILED ', NFAIL_REZ, 'time(s)'
  657. WRITE(*,*) 'Max residual computing test adjusted error measure was ', TMP_REZ
  658. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  659. NFAIL_TOTAL = NFAIL_TOTAL + NFAIL_REZ
  660. END IF
  661. IF ( NFAIL_TOTAL == 0 ) THEN
  662. WRITE(*,*) '>>>> SGEDMD :: ALL TESTS PASSED.'
  663. ELSE
  664. WRITE(*,*) NFAIL_TOTAL, 'FAILURES!'
  665. WRITE(*,*) '>>>>>>>>>>>>>> SGEDMD :: TESTS FAILED. CHECK THE IMPLEMENTATION.'
  666. END IF
  667. IF ( TEST_QRDMD ) THEN
  668. WRITE(*,*)
  669. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  670. WRITE(*,*) ' Test summary for SGEDMDQ :'
  671. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  672. WRITE(*,*)
  673. IF ( NFAIL_SVDIFF == 0 ) THEN
  674. WRITE(*,*) '>>>> SGEDMD and SGEDMDQ computed singular &
  675. &values test PASSED.'
  676. ELSE
  677. WRITE(*,*) 'SGEDMD and SGEDMDQ discrepancies in &
  678. &the singular values unacceptable ', &
  679. NFAIL_SVDIFF, ' times. Test FAILED.'
  680. WRITE(*,*) 'The maximal discrepancy in the singular values (relative to the norm) was ', SVDIFF
  681. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  682. NFAILQ_TOTAL = NFAILQ_TOTAL + NFAIL_SVDIFF
  683. END IF
  684. IF ( NFAIL_F_QR == 0 ) THEN
  685. WRITE(*,*) '>>>> F - Q*R test PASSED.'
  686. ELSE
  687. WRITE(*,*) 'F - Q*R test FAILED ', NFAIL_F_QR, ' time(s)'
  688. WRITE(*,*) 'The largest relative residual was ', TMP_FQR
  689. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  690. NFAILQ_TOTAL = NFAILQ_TOTAL + NFAIL_F_QR
  691. END IF
  692. IF ( NFAIL_REZQ == 0 ) THEN
  693. WRITE(*,*) '>>>> Rezidual computation test PASSED.'
  694. ELSE
  695. WRITE(*,*) 'Rezidual computation test FAILED ', NFAIL_REZQ, 'time(s)'
  696. WRITE(*,*) 'Max residual computing test adjusted error measure was ', TMP_REZQ
  697. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  698. NFAILQ_TOTAL = NFAILQ_TOTAL + NFAIL_REZQ
  699. END IF
  700. IF ( NFAILQ_TOTAL == 0 ) THEN
  701. WRITE(*,*) '>>>>>>> SGEDMDQ :: ALL TESTS PASSED.'
  702. ELSE
  703. WRITE(*,*) NFAILQ_TOTAL, 'FAILURES!'
  704. WRITE(*,*) '>>>>>>> SGEDMDQ :: TESTS FAILED. CHECK THE IMPLEMENTATION.'
  705. END IF
  706. END IF
  707. WRITE(*,*)
  708. WRITE(*,*) 'Test completed.'
  709. STOP
  710. END