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.

dchkdmd.f90 30 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813
  1. ! This is a test program for checking the implementations of
  2. ! the implementations of the following subroutines
  3. !
  4. ! DGEDMD for computation of the
  5. ! Dynamic Mode Decomposition (DMD)
  6. ! DGEDMDQ 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: real64
  63. IMPLICIT NONE
  64. integer, parameter :: WP = real64
  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, j, 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 DAXPY, DGEEV, DGEMM, DGEMV, DLACPY, DLASCL
  96. EXTERNAL DLARNV, DLATMR
  97. !.....external subroutines DMD package, part 1
  98. ! subroutines under test
  99. EXTERNAL DGEDMD, DGEDMDQ
  100. !..... external functions (BLAS and LAPACK)
  101. EXTERNAL DLAMCH, DLANGE, DNRM2
  102. REAL(KIND=WP) :: DLAMCH, DLANGE, DNRM2
  103. EXTERNAL LSAME
  104. LOGICAL LSAME
  105. INTRINSIC ABS, INT, MIN, MAX
  106. !............................................................
  107. ! The test is always in pairs : ( DGEDMD and DGEDMDQ )
  108. ! because the test includes comparing the results (in pairs).
  109. !.....................................................................................
  110. TEST_QRDMD = .TRUE. ! This code by default performs tests on DGEDMDQ
  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 = DLAMCH( 'P' ) ! machine precision DP
  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 DLATMR( 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 DGEEV( '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 DLASCL( 'G', 0, 0, TMP, ONE, M, M, A, M, INFO )
  242. CALL DLASCL( 'G', 0, 0, TMP, ONE, M, 2, EIGA, M, INFO )
  243. ! Compute the norm of A
  244. ANORM = DLANGE( 'F', N, N, A, M, WDUMMY )
  245. IF ( K_TRAJ == 2 ) THEN
  246. ! generate data with two inital conditions
  247. CALL DLARNV(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 DGEMV( '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 DLARNV(2, ISEED, M, F1(1,1) )
  256. DO i = 1, N-N/2
  257. CALL DGEMV( '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. CALL DLARNV(2, ISEED, M, F(1,1) )
  264. DO i = 1, N
  265. CALL DGEMV( 'N', M, M, ONE, A, M, F(1,i), 1, ZERO, &
  266. F(1,i+1), 1 )
  267. END DO
  268. X0(1:M,1:N) = F(1:M,1:N)
  269. Y0(1:M,1:N) = F(1:M,2:N+1)
  270. END IF
  271. XNORM = DLANGE( 'F', M, N, X0, LDX, WDUMMY )
  272. YNORM = DLANGE( 'F', M, N, Y0, LDX, WDUMMY )
  273. !............................................................
  274. DO iJOBZ = 1, 4
  275. SELECT CASE ( iJOBZ )
  276. CASE(1)
  277. JOBZ = 'V' ! Ritz vectors will be computed
  278. RESIDS = 'R' ! Residuals will be computed
  279. CASE(2)
  280. JOBZ = 'V'
  281. RESIDS = 'N'
  282. CASE(3)
  283. JOBZ = 'F' ! Ritz vectors in factored form
  284. RESIDS = 'N'
  285. CASE(4)
  286. JOBZ = 'N'
  287. RESIDS = 'N'
  288. END SELECT
  289. DO iJOBREF = 1, 3
  290. SELECT CASE ( iJOBREF )
  291. CASE(1)
  292. JOBREF = 'R' ! Data for refined Ritz vectors
  293. CASE(2)
  294. JOBREF = 'E' ! Exact DMD vectors
  295. CASE(3)
  296. JOBREF = 'N'
  297. END SELECT
  298. DO iSCALE = 1, 4
  299. SELECT CASE ( iSCALE )
  300. CASE(1)
  301. SCALE = 'S' ! X data normalized
  302. CASE(2)
  303. SCALE = 'C' ! X normalized, consist. check
  304. CASE(3)
  305. SCALE = 'Y' ! Y data normalized
  306. CASE(4)
  307. SCALE = 'N'
  308. END SELECT
  309. DO iNRNK = -1, -2, -1
  310. ! Two truncation strategies. The "-2" case for R&D
  311. ! purposes only - it uses possibly low accuracy small
  312. ! singular values, in which case the formulas used in
  313. ! the DMD are highly sensitive.
  314. NRNK = iNRNK
  315. DO iWHTSVD = 1, 4
  316. ! Check all four options to compute the POD basis
  317. ! via the SVD.
  318. WHTSVD = iWHTSVD
  319. DO LWMINOPT = 1, 2
  320. ! Workspace query for the minimal (1) and for the optimal
  321. ! (2) workspace lengths determined by workspace query.
  322. X(1:M,1:N) = X0(1:M,1:N)
  323. Y(1:M,1:N) = Y0(1:M,1:N)
  324. ! DGEDMD: Workspace query and workspace allocation
  325. CALL DGEDMD( SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, M, &
  326. N, X, LDX, Y, LDY, NRNK, TOL, K, REIG, IEIG, Z, &
  327. LDZ, RES, AU, LDAU, W, LDW, S, LDS, WDUMMY, -1, &
  328. IDUMMY, -1, INFO )
  329. LIWORK = IDUMMY(1)
  330. ALLOCATE( IWORK(LIWORK) )
  331. LWORK = INT(WDUMMY(LWMINOPT))
  332. ALLOCATE( WORK(LWORK) )
  333. ! DGEDMD test: CALL DGEDMD
  334. CALL DGEDMD( SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, M, &
  335. N, X, LDX, Y, LDY, NRNK, TOL, K, REIG, IEIG, Z, &
  336. LDZ, RES, AU, LDAU, W, LDW, S, LDS, WORK, LWORK,&
  337. IWORK, LIWORK, INFO )
  338. SINGVX(1:N) = WORK(1:N)
  339. !...... DGEDMD check point
  340. IF ( LSAME(JOBZ,'V') ) THEN
  341. ! Check that Z = X*W, on return from DGEDMD
  342. ! This checks that the returned aigenvectors in Z are
  343. ! the product of the SVD'POD basis returned in X
  344. ! and the eigenvectors of the rayleigh quotient
  345. ! returned in W
  346. CALL DGEMM( 'N', 'N', M, K, K, ONE, X, LDX, W, LDW, &
  347. ZERO, Z1, LDZ )
  348. TMP = ZERO
  349. DO i = 1, K
  350. CALL DAXPY( M, -ONE, Z(1,i), 1, Z1(1,i), 1)
  351. TMP = MAX(TMP, DNRM2( M, Z1(1,i), 1 ) )
  352. END DO
  353. TMP_ZXW = MAX(TMP_ZXW, TMP )
  354. IF ( TMP_ZXW > 10*M*EPS ) THEN
  355. NFAIL_Z_XV = NFAIL_Z_XV + 1
  356. WRITE(*,*) ':( .................DGEDMD FAILED!', &
  357. 'Check the code for implementation errors.'
  358. WRITE(*,*) 'The input parameters were ',&
  359. SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, &
  360. M, N, LDX, LDY, NRNK, TOL
  361. END IF
  362. END IF
  363. !...... DGEDMD check point
  364. IF ( LSAME(JOBREF,'R') ) THEN
  365. ! The matrix A*U is returned for computing refined Ritz vectors.
  366. ! Check that A*U is computed correctly using the formula
  367. ! A*U = Y * V * inv(SIGMA). This depends on the
  368. ! accuracy in the computed singular values and vectors of X.
  369. ! See the paper for an error analysis.
  370. ! Note that the left singular vectors of the input matrix X
  371. ! are returned in the array X.
  372. CALL DGEMM( 'N', 'N', M, K, M, ONE, A, LDA, X, LDX, &
  373. ZERO, Z1, LDZ )
  374. TMP = ZERO
  375. DO i = 1, K
  376. CALL DAXPY( M, -ONE, AU(1,i), 1, Z1(1,i), 1)
  377. TMP = MAX( TMP, DNRM2( M, Z1(1,i),1 ) * &
  378. SINGVX(K)/(ANORM*SINGVX(1)) )
  379. END DO
  380. TMP_AU = MAX( TMP_AU, TMP )
  381. IF ( TMP > TOL2 ) THEN
  382. NFAIL_AU = NFAIL_AU + 1
  383. WRITE(*,*) ':( .................DGEDMD FAILED!', &
  384. 'Check the code for implementation errors.'
  385. WRITE(*,*) 'The input parameters were ',&
  386. SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, &
  387. M, N, LDX, LDY, NRNK, TOL
  388. END IF
  389. ELSEIF ( LSAME(JOBREF,'E') ) THEN
  390. ! The unscaled vectors of the Exact DMD are computed.
  391. ! This option is included for the sake of completeness,
  392. ! for users who prefer the Exact DMD vectors. The
  393. ! returned vectors are in the real form, in the same way
  394. ! as the Ritz vectors. Here we just save the vectors
  395. ! and test them separately using a Matlab script.
  396. CALL DGEMM( 'N', 'N', M, K, M, ONE, A, LDA, AU, LDAU, ZERO, Y1, M )
  397. i=1
  398. DO WHILE ( i <= K )
  399. IF ( IEIG(i) == ZERO ) THEN
  400. ! have a real eigenvalue with real eigenvector
  401. CALL DAXPY( M, -REIG(i), AU(1,i), 1, Y1(1,i), 1 )
  402. RESEX(i) = DNRM2( M, Y1(1,i), 1) / DNRM2(M,AU(1,i),1)
  403. i = i + 1
  404. ELSE
  405. ! Have a complex conjugate pair
  406. ! REIG(i) +- sqrt(-1)*IMEIG(i).
  407. ! Since all computation is done in real
  408. ! arithmetic, the formula for the residual
  409. ! is recast for real representation of the
  410. ! complex conjugate eigenpair. See the
  411. ! description of RES.
  412. AB(1,1) = REIG(i)
  413. AB(2,1) = -IEIG(i)
  414. AB(1,2) = IEIG(i)
  415. AB(2,2) = REIG(i)
  416. CALL DGEMM( 'N', 'N', M, 2, 2, -ONE, AU(1,i), &
  417. M, AB, 2, ONE, Y1(1,i), M )
  418. RESEX(i) = DLANGE( 'F', M, 2, Y1(1,i), M, &
  419. WORK )/ DLANGE( 'F', M, 2, AU(1,i), M, &
  420. WORK )
  421. RESEX(i+1) = RESEX(i)
  422. i = i + 2
  423. END IF
  424. END DO
  425. END IF
  426. !...... DGEDMD check point
  427. IF ( LSAME(RESIDS, 'R') ) THEN
  428. ! Compare the residuals returned by DGEDMD with the
  429. ! explicitly computed residuals using the matrix A.
  430. ! Compute explicitly Y1 = A*Z
  431. CALL DGEMM( 'N', 'N', M, K, M, ONE, A, LDA, Z, LDZ, ZERO, Y1, M )
  432. ! ... and then A*Z(:,i) - LAMBDA(i)*Z(:,i), using the real forms
  433. ! of the invariant subspaces that correspond to complex conjugate
  434. ! pairs of eigencalues. (See the description of Z in DGEDMD,)
  435. i = 1
  436. DO WHILE ( i <= K )
  437. IF ( IEIG(i) == ZERO ) THEN
  438. ! have a real eigenvalue with real eigenvector
  439. CALL DAXPY( M, -REIG(i), Z(1,i), 1, Y1(1,i), 1 )
  440. RES1(i) = DNRM2( M, Y1(1,i), 1)
  441. i = i + 1
  442. ELSE
  443. ! Have a complex conjugate pair
  444. ! REIG(i) +- sqrt(-1)*IMEIG(i).
  445. ! Since all computation is done in real
  446. ! arithmetic, the formula for the residual
  447. ! is recast for real representation of the
  448. ! complex conjugate eigenpair. See the
  449. ! description of RES.
  450. AB(1,1) = REIG(i)
  451. AB(2,1) = -IEIG(i)
  452. AB(1,2) = IEIG(i)
  453. AB(2,2) = REIG(i)
  454. CALL DGEMM( 'N', 'N', M, 2, 2, -ONE, Z(1,i), &
  455. M, AB, 2, ONE, Y1(1,i), M )
  456. RES1(i) = DLANGE( 'F', M, 2, Y1(1,i), M, &
  457. WORK )
  458. RES1(i+1) = RES1(i)
  459. i = i + 2
  460. END IF
  461. END DO
  462. TMP = ZERO
  463. DO i = 1, K
  464. TMP = MAX( TMP, ABS(RES(i) - RES1(i)) * &
  465. SINGVX(K)/(ANORM*SINGVX(1)) )
  466. END DO
  467. TMP_REZ = MAX( TMP_REZ, TMP )
  468. IF ( TMP > TOL2 ) THEN
  469. NFAIL_REZ = NFAIL_REZ + 1
  470. WRITE(*,*) ':( ..................DGEDMD FAILED!', &
  471. 'Check the code for implementation errors.'
  472. WRITE(*,*) 'The input parameters were ',&
  473. SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, &
  474. M, N, LDX, LDY, NRNK, TOL
  475. END IF
  476. IF ( LSAME(JOBREF,'E') ) THEN
  477. TMP = ZERO
  478. DO i = 1, K
  479. TMP = MAX( TMP, ABS(RES1(i) - RESEX(i))/(RES1(i)+RESEX(i)) )
  480. END DO
  481. TMP_EX = MAX(TMP_EX,TMP)
  482. END IF
  483. END IF
  484. !..... store the results for inspection
  485. DO i = 1, K
  486. LAMBDA(i,1) = REIG(i)
  487. LAMBDA(i,2) = IEIG(i)
  488. END DO
  489. DEALLOCATE(IWORK)
  490. DEALLOCATE(WORK)
  491. !======================================================================
  492. ! Now test the DGEDMDQ
  493. !======================================================================
  494. IF ( TEST_QRDMD .AND. (K_TRAJ == 1) ) THEN
  495. RJOBDATA(2) = 1
  496. F1 = F
  497. ! DGEDMDQ test: Workspace query and workspace allocation
  498. CALL DGEDMDQ( SCALE, JOBZ, RESIDS, WANTQ, WANTR, &
  499. JOBREF, WHTSVD, M, N+1, F1, LDF, X, LDX, Y, &
  500. LDY, NRNK, TOL, KQ, REIGQ, IEIGQ, Z, LDZ, &
  501. RES, AU, LDAU, W, LDW, S, LDS, WDUMMY, &
  502. -1, IDUMMY, -1, INFO )
  503. LIWORK = IDUMMY(1)
  504. ALLOCATE( IWORK(LIWORK) )
  505. LWORK = INT(WDUMMY(LWMINOPT))
  506. ALLOCATE(WORK(LWORK))
  507. ! DGEDMDQ test: CALL DGEDMDQ
  508. CALL DGEDMDQ( SCALE, JOBZ, RESIDS, WANTQ, WANTR, &
  509. JOBREF, WHTSVD, M, N+1, F1, LDF, X, LDX, Y, &
  510. LDY, NRNK, TOL, KQ, REIGQ, IEIGQ, Z, LDZ, &
  511. RES, AU, LDAU, W, LDW, S, LDS, &
  512. WORK, LWORK, IWORK, LIWORK, INFO )
  513. SINGVQX(1:KQ) = WORK(MIN(M,N+1)+1: MIN(M,N+1)+KQ)
  514. !..... DGEDMDQ check point
  515. IF ( KQ /= K ) THEN
  516. KDIFF = KDIFF+1
  517. END IF
  518. TMP = ZERO
  519. DO i = 1, MIN(K, KQ)
  520. TMP = MAX(TMP, ABS(SINGVX(i)-SINGVQX(i)) / &
  521. SINGVX(1) )
  522. END DO
  523. SVDIFF = MAX( SVDIFF, TMP )
  524. IF ( TMP > M*N*EPS ) THEN
  525. WRITE(*,*) 'FAILED! Something was wrong with the run.'
  526. NFAIL_SVDIFF = NFAIL_SVDIFF + 1
  527. DO j =1, 3
  528. write(*,*) j, SINGVX(j), SINGVQX(j)
  529. read(*,*)
  530. END DO
  531. END IF
  532. !..... DGEDMDQ check point
  533. IF ( LSAME(WANTQ,'Q') .AND. LSAME(WANTR,'R') ) THEN
  534. ! Check that the QR factors are computed and returned
  535. ! as requested. The residual ||F-Q*R||_F / ||F||_F
  536. ! is compared to M*N*EPS.
  537. F2 = F
  538. CALL DGEMM( 'N', 'N', M, N+1, MIN(M,N+1), -ONE, F1, &
  539. LDF, Y, LDY, ONE, F2, LDF )
  540. TMP_FQR = DLANGE( 'F', M, N+1, F2, LDF, WORK ) / &
  541. DLANGE( 'F', M, N+1, F, LDF, WORK )
  542. IF ( TMP_FQR > TOL2 ) THEN
  543. WRITE(*,*) 'FAILED! Something was wrong with the run.'
  544. NFAIL_F_QR = NFAIL_F_QR + 1
  545. END IF
  546. END IF
  547. !..... DGEDMDQ check point
  548. IF ( LSAME(RESIDS, 'R') ) THEN
  549. ! Compare the residuals returned by DGEDMDQ with the
  550. ! explicitly computed residuals using the matrix A.
  551. ! Compute explicitly Y1 = A*Z
  552. CALL DGEMM( 'N', 'N', M, KQ, M, ONE, A, M, Z, M, ZERO, Y1, M )
  553. ! ... and then A*Z(:,i) - LAMBDA(i)*Z(:,i), using the real forms
  554. ! of the invariant subspaces that correspond to complex conjugate
  555. ! pairs of eigencalues. (See the description of Z in DGEDMDQ)
  556. i = 1
  557. DO WHILE ( i <= KQ )
  558. IF ( IEIGQ(i) == ZERO ) THEN
  559. ! have a real eigenvalue with real eigenvector
  560. CALL DAXPY( M, -REIGQ(i), Z(1,i), 1, Y1(1,i), 1 )
  561. ! Y(1:M,i) = Y(1:M,i) - REIG(i)*Z(1:M,i)
  562. RES1(i) = DNRM2( M, Y1(1,i), 1)
  563. i = i + 1
  564. ELSE
  565. ! Have a complex conjugate pair
  566. ! REIG(i) +- sqrt(-1)*IMEIG(i).
  567. ! Since all computation is done in real
  568. ! arithmetic, the formula for the residual
  569. ! is recast for real representation of the
  570. ! complex conjugate eigenpair. See the
  571. ! description of RES.
  572. AB(1,1) = REIGQ(i)
  573. AB(2,1) = -IEIGQ(i)
  574. AB(1,2) = IEIGQ(i)
  575. AB(2,2) = REIGQ(i)
  576. CALL DGEMM( 'N', 'N', M, 2, 2, -ONE, Z(1,i), &
  577. M, AB, 2, ONE, Y1(1,i), M ) ! BLAS CALL
  578. ! Y(1:M,i:i+1) = Y(1:M,i:i+1) - Z(1:M,i:i+1) * AB ! INTRINSIC
  579. RES1(i) = DLANGE( 'F', M, 2, Y1(1,i), M, &
  580. WORK ) ! LAPACK CALL
  581. RES1(i+1) = RES1(i)
  582. i = i + 2
  583. END IF
  584. END DO
  585. TMP = ZERO
  586. DO i = 1, KQ
  587. TMP = MAX( TMP, ABS(RES(i) - RES1(i)) * &
  588. SINGVQX(K)/(ANORM*SINGVQX(1)) )
  589. END DO
  590. TMP_REZQ = MAX( TMP_REZQ, TMP )
  591. IF ( TMP > TOL2 ) THEN
  592. NFAIL_REZQ = NFAIL_REZQ + 1
  593. WRITE(*,*) '................ DGEDMDQ FAILED!', &
  594. 'Check the code for implementation errors.'
  595. STOP
  596. END IF
  597. END IF
  598. DO i = 1, KQ
  599. LAMBDAQ(i,1) = REIGQ(i)
  600. LAMBDAQ(i,2) = IEIGQ(i)
  601. END DO
  602. DEALLOCATE(WORK)
  603. DEALLOCATE(IWORK)
  604. END IF ! TEST_QRDMD
  605. !======================================================================
  606. END DO ! LWMINOPT
  607. !write(*,*) 'LWMINOPT loop completed'
  608. END DO ! WHTSVD LOOP
  609. !write(*,*) 'WHTSVD loop completed'
  610. END DO ! NRNK LOOP
  611. !write(*,*) 'NRNK loop completed'
  612. END DO ! SCALE LOOP
  613. !write(*,*) 'SCALE loop completed'
  614. END DO ! JOBF LOOP
  615. !write(*,*) 'JOBREF loop completed'
  616. END DO ! JOBZ LOOP
  617. !write(*,*) 'JOBZ loop completed'
  618. END DO ! MODE -6:6
  619. !write(*,*) 'MODE loop completed'
  620. END DO ! 1 or 2 trajectories
  621. !write(*,*) 'trajectories loop completed'
  622. DEALLOCATE(A)
  623. DEALLOCATE(AC)
  624. DEALLOCATE(DA)
  625. DEALLOCATE(DL)
  626. DEALLOCATE(F)
  627. DEALLOCATE(F1)
  628. DEALLOCATE(F2)
  629. DEALLOCATE(X)
  630. DEALLOCATE(X0)
  631. DEALLOCATE(SINGVX)
  632. DEALLOCATE(SINGVQX)
  633. DEALLOCATE(Y)
  634. DEALLOCATE(Y0)
  635. DEALLOCATE(Y1)
  636. DEALLOCATE(Z)
  637. DEALLOCATE(Z1)
  638. DEALLOCATE(RES)
  639. DEALLOCATE(RES1)
  640. DEALLOCATE(RESEX)
  641. DEALLOCATE(REIG)
  642. DEALLOCATE(IEIG)
  643. DEALLOCATE(REIGQ)
  644. DEALLOCATE(IEIGQ)
  645. DEALLOCATE(REIGA)
  646. DEALLOCATE(IEIGA)
  647. DEALLOCATE(VA)
  648. DEALLOCATE(LAMBDA)
  649. DEALLOCATE(LAMBDAQ)
  650. DEALLOCATE(EIGA)
  651. DEALLOCATE(W)
  652. DEALLOCATE(AU)
  653. DEALLOCATE(S)
  654. !............................................................
  655. ! Generate random M-by-M matrix A. Use DLATMR from
  656. END DO ! LLOOP
  657. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  658. WRITE(*,*) ' Test summary for DGEDMD :'
  659. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  660. WRITE(*,*)
  661. IF ( NFAIL_Z_XV == 0 ) THEN
  662. WRITE(*,*) '>>>> Z - U*V test PASSED.'
  663. ELSE
  664. WRITE(*,*) 'Z - U*V test FAILED ', NFAIL_Z_XV, ' time(s)'
  665. WRITE(*,*) 'Max error ||Z-U*V||_F was ', TMP_ZXW
  666. NFAIL_TOTAL = NFAIL_TOTAL + NFAIL_Z_XV
  667. END IF
  668. IF ( NFAIL_AU == 0 ) THEN
  669. WRITE(*,*) '>>>> A*U test PASSED. '
  670. ELSE
  671. WRITE(*,*) 'A*U test FAILED ', NFAIL_AU, ' time(s)'
  672. WRITE(*,*) 'Max A*U test adjusted error measure was ', TMP_AU
  673. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  674. NFAIL_TOTAL = NFAIL_TOTAL + NFAIL_AU
  675. END IF
  676. IF ( NFAIL_REZ == 0 ) THEN
  677. WRITE(*,*) '>>>> Rezidual computation test PASSED.'
  678. ELSE
  679. WRITE(*,*) 'Rezidual computation test FAILED ', NFAIL_REZ, 'time(s)'
  680. WRITE(*,*) 'Max residual computing test adjusted error measure was ', TMP_REZ
  681. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  682. NFAIL_TOTAL = NFAIL_TOTAL + NFAIL_REZ
  683. END IF
  684. IF ( NFAIL_TOTAL == 0 ) THEN
  685. WRITE(*,*) '>>>> DGEDMD :: ALL TESTS PASSED.'
  686. ELSE
  687. WRITE(*,*) NFAIL_TOTAL, 'FAILURES!'
  688. WRITE(*,*) '>>>>>>>>>>>>>> DGEDMD :: TESTS FAILED. CHECK THE IMPLEMENTATION.'
  689. END IF
  690. IF ( TEST_QRDMD ) THEN
  691. WRITE(*,*)
  692. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  693. WRITE(*,*) ' Test summary for DGEDMDQ :'
  694. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  695. WRITE(*,*)
  696. IF ( NFAIL_SVDIFF == 0 ) THEN
  697. WRITE(*,*) '>>>> DGEDMD and DGEDMDQ computed singular &
  698. &values test PASSED.'
  699. ELSE
  700. WRITE(*,*) 'DGEDMD and DGEDMDQ discrepancies in &
  701. &the singular values unacceptable ', &
  702. NFAIL_SVDIFF, ' times. Test FAILED.'
  703. WRITE(*,*) 'The maximal discrepancy in the singular values (relative to the norm) was ', SVDIFF
  704. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  705. NFAILQ_TOTAL = NFAILQ_TOTAL + NFAIL_SVDIFF
  706. END IF
  707. IF ( NFAIL_F_QR == 0 ) THEN
  708. WRITE(*,*) '>>>> F - Q*R test PASSED.'
  709. ELSE
  710. WRITE(*,*) 'F - Q*R test FAILED ', NFAIL_F_QR, ' time(s)'
  711. WRITE(*,*) 'The largest relative residual was ', TMP_FQR
  712. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  713. NFAILQ_TOTAL = NFAILQ_TOTAL + NFAIL_F_QR
  714. END IF
  715. IF ( NFAIL_REZQ == 0 ) THEN
  716. WRITE(*,*) '>>>> Rezidual computation test PASSED.'
  717. ELSE
  718. WRITE(*,*) 'Rezidual computation test FAILED ', NFAIL_REZQ, 'time(s)'
  719. WRITE(*,*) 'Max residual computing test adjusted error measure was ', TMP_REZQ
  720. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  721. NFAILQ_TOTAL = NFAILQ_TOTAL + NFAIL_REZQ
  722. END IF
  723. IF ( NFAILQ_TOTAL == 0 ) THEN
  724. WRITE(*,*) '>>>>>>> DGEDMDQ :: ALL TESTS PASSED.'
  725. ELSE
  726. WRITE(*,*) NFAILQ_TOTAL, 'FAILURES!'
  727. WRITE(*,*) '>>>>>>> DGEDMDQ :: TESTS FAILED. CHECK THE IMPLEMENTATION.'
  728. END IF
  729. END IF
  730. WRITE(*,*)
  731. WRITE(*,*) 'Test completed.'
  732. STOP
  733. END