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.

zchkdmd.f90 27 kB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745
  1. ! This is a test program for checking the implementations of
  2. ! the implementations of the following subroutines
  3. !
  4. ! ZGEDMD, for computation of the
  5. ! Dynamic Mode Decomposition (DMD)
  6. ! ZGEDMDQ, 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. !............................................................
  39. !............................................................
  40. !
  41. PROGRAM DMD_TEST
  42. use iso_fortran_env, only: real64
  43. IMPLICIT NONE
  44. integer, parameter :: WP = real64
  45. !............................................................
  46. REAL(KIND=WP), PARAMETER :: ONE = 1.0_WP
  47. REAL(KIND=WP), PARAMETER :: ZERO = 0.0_WP
  48. COMPLEX(KIND=WP), PARAMETER :: ZONE = ( 1.0_WP, 0.0_WP )
  49. COMPLEX(KIND=WP), PARAMETER :: ZZERO = ( 0.0_WP, 0.0_WP )
  50. !............................................................
  51. REAL(KIND=WP), ALLOCATABLE, DIMENSION(:) :: RES, &
  52. RES1, RESEX, SINGVX, SINGVQX, WORK
  53. INTEGER , ALLOCATABLE, DIMENSION(:) :: IWORK
  54. REAL(KIND=WP) :: WDUMMY(2)
  55. INTEGER :: IDUMMY(4), ISEED(4)
  56. REAL(KIND=WP) :: ANORM, COND, CONDL, CONDR, EPS, &
  57. TOL, TOL2, SVDIFF, TMP, TMP_AU, &
  58. TMP_FQR, TMP_REZ, TMP_REZQ, TMP_ZXW, &
  59. TMP_EX
  60. !............................................................
  61. COMPLEX(KIND=WP) :: ZMAX
  62. INTEGER :: LZWORK
  63. COMPLEX(KIND=WP), ALLOCATABLE, DIMENSION(:,:) :: ZA, ZAC, &
  64. ZAU, ZF, ZF0, ZF1, ZS, ZW, &
  65. ZX, ZX0, ZY, ZY0, ZY1, ZZ, ZZ1
  66. COMPLEX(KIND=WP), ALLOCATABLE, DIMENSION(:) :: ZDA, ZDR, &
  67. ZDL, ZEIGS, ZEIGSA, ZWORK
  68. COMPLEX(KIND=WP) :: ZDUMMY(22), ZDUM2X2(2,2)
  69. !............................................................
  70. INTEGER :: K, KQ, LDF, LDS, LDA, LDAU, LDW, LDX, LDY, &
  71. LDZ, LIWORK, LWORK, M, N, LLOOP, NRNK, NRNKsp
  72. INTEGER :: i, iJOBREF, iJOBZ, iSCALE, INFO, j, &
  73. NFAIL, NFAIL_AU, NFAIL_F_QR, NFAIL_REZ, &
  74. NFAIL_REZQ, NFAIL_SVDIFF, NFAIL_TOTAL, NFAILQ_TOTAL, &
  75. NFAIL_Z_XV, MODE, MODEL, MODER, WHTSVD, &
  76. WHTSVDsp
  77. INTEGER :: iNRNK, iWHTSVD, K_TRAJ, LWMINOPT
  78. CHARACTER :: GRADE, JOBREF, JOBZ, PIVTNG, RSIGN, &
  79. SCALE, RESIDS, WANTQ, WANTR
  80. LOGICAL :: TEST_QRDMD
  81. !.....external subroutines (BLAS and LAPACK)
  82. EXTERNAL DAXPY, DGEEV, DGEMM, DGEMV, DLACPY, DLASCL
  83. EXTERNAL ZGEEV, ZGEMV, ZLASCL
  84. EXTERNAL ZLARNV, ZLATMR
  85. EXTERNAL ZAXPY, ZGEMM
  86. !.....external subroutines DMD package, part 1
  87. ! subroutines under test
  88. EXTERNAL ZGEDMD, ZGEDMDQ
  89. !.....external functions (BLAS and LAPACK)
  90. EXTERNAL DLAMCH, DZNRM2
  91. REAL(KIND=WP) :: DLAMCH, DZNRM2
  92. REAL(KIND=WP) :: ZLANGE
  93. EXTERNAL IZAMAX
  94. INTEGER IZAMAX
  95. EXTERNAL LSAME
  96. LOGICAL LSAME
  97. INTRINSIC ABS, INT, MIN, MAX, SIGN
  98. !............................................................
  99. ! The test is always in pairs : ( ZGEDMD and ZGEDMDQ )
  100. ! because the test includes comparing the results (in pairs).
  101. !.....................................................................................
  102. TEST_QRDMD = .TRUE. ! This code by default performs tests on ZGEDMDQ
  103. ! Since the QR factorizations based algorithm is designed for
  104. ! single trajectory data, only single trajectory tests will
  105. ! be performed with xGEDMDQ.
  106. WANTQ = 'Q'
  107. WANTR = 'R'
  108. !.................................................................................
  109. EPS = DLAMCH( 'P' ) ! machine precision DP
  110. ! Global counters of failures of some particular tests
  111. NFAIL = 0
  112. NFAIL_REZ = 0
  113. NFAIL_REZQ = 0
  114. NFAIL_Z_XV = 0
  115. NFAIL_F_QR = 0
  116. NFAIL_AU = 0
  117. NFAIL_SVDIFF = 0
  118. NFAIL_TOTAL = 0
  119. NFAILQ_TOTAL = 0
  120. DO LLOOP = 1, 4
  121. WRITE(*,*) 'L Loop Index = ', LLOOP
  122. ! Set the dimensions of the problem ...
  123. WRITE(*,*) 'M = '
  124. READ(*,*) M
  125. WRITE(*,*) M
  126. ! ... and the number of snapshots.
  127. WRITE(*,*) 'N = '
  128. READ(*,*) N
  129. WRITE(*,*) N
  130. ! ... Test the dimensions
  131. IF ( ( MIN(M,N) == 0 ) .OR. ( M < N ) ) THEN
  132. WRITE(*,*) 'Bad dimensions. Required: M >= N > 0.'
  133. STOP
  134. END IF
  135. !.............
  136. ! The seed inside the LLOOP so that each pass can be reproduced easily.
  137. ISEED(1) = 4
  138. ISEED(2) = 3
  139. ISEED(3) = 2
  140. ISEED(4) = 1
  141. LDA = M
  142. LDF = M
  143. LDX = M
  144. LDY = M
  145. LDW = N
  146. LDZ = M
  147. LDAU = M
  148. LDS = N
  149. TMP_ZXW = ZERO
  150. TMP_AU = ZERO
  151. TMP_REZ = ZERO
  152. TMP_REZQ = ZERO
  153. SVDIFF = ZERO
  154. TMP_EX = ZERO
  155. ALLOCATE( ZA(LDA,M) )
  156. ALLOCATE( ZAC(LDA,M) )
  157. ALLOCATE( ZF(LDF,N+1) )
  158. ALLOCATE( ZF0(LDF,N+1) )
  159. ALLOCATE( ZF1(LDF,N+1) )
  160. ALLOCATE( ZX(LDX,N) )
  161. ALLOCATE( ZX0(LDX,N) )
  162. ALLOCATE( ZY(LDY,N+1) )
  163. ALLOCATE( ZY0(LDY,N+1) )
  164. ALLOCATE( ZY1(LDY,N+1) )
  165. ALLOCATE( ZAU(LDAU,N) )
  166. ALLOCATE( ZW(LDW,N) )
  167. ALLOCATE( ZS(LDS,N) )
  168. ALLOCATE( ZZ(LDZ,N) )
  169. ALLOCATE( ZZ1(LDZ,N) )
  170. ALLOCATE( RES(N) )
  171. ALLOCATE( RES1(N) )
  172. ALLOCATE( RESEX(N) )
  173. ALLOCATE( ZEIGS(N) )
  174. ALLOCATE( SINGVX(N) )
  175. ALLOCATE( SINGVQX(N) )
  176. TOL = M*EPS
  177. ! This mimics O(M*N)*EPS bound for accumulated roundoff error.
  178. ! The factor 10 is somewhat arbitrary.
  179. TOL2 = 10*M*N*EPS
  180. !.............
  181. DO K_TRAJ = 1, 2
  182. ! Number of intial conditions in the simulation/trajectories (1 or 2)
  183. COND = 1.0D4
  184. ZMAX = (1.0D1,1.0D1)
  185. RSIGN = 'F'
  186. GRADE = 'N'
  187. MODEL = 6
  188. CONDL = 1.0D1
  189. MODER = 6
  190. CONDR = 1.0D1
  191. PIVTNG = 'N'
  192. ! Loop over all parameter MODE values for ZLATMR (+1,..,+6)
  193. DO MODE = 1, 6
  194. ALLOCATE( IWORK(2*M) )
  195. ALLOCATE( ZDA(M) )
  196. ALLOCATE( ZDL(M) )
  197. ALLOCATE( ZDR(M) )
  198. CALL ZLATMR( M, M, 'N', ISEED, 'N', ZDA, MODE, COND, &
  199. ZMAX, RSIGN, GRADE, ZDL, MODEL, CONDL, &
  200. ZDR, MODER, CONDR, PIVTNG, IWORK, M, M, &
  201. ZERO, -ONE, 'N', ZA, LDA, IWORK(M+1), INFO )
  202. DEALLOCATE( ZDR )
  203. DEALLOCATE( ZDL )
  204. DEALLOCATE( ZDA )
  205. DEALLOCATE( IWORK )
  206. LZWORK = MAX(1,2*M)
  207. ALLOCATE( ZEIGSA(M) )
  208. ALLOCATE( ZWORK(LZWORK) )
  209. ALLOCATE( WORK(2*M) )
  210. ZAC(1:M,1:M) = ZA(1:M,1:M)
  211. CALL ZGEEV( 'N','N', M, ZAC, LDA, ZEIGSA, ZDUM2X2, 2, &
  212. ZDUM2X2, 2, ZWORK, LZWORK, WORK, INFO ) ! LAPACK CALL
  213. DEALLOCATE(WORK)
  214. DEALLOCATE(ZWORK)
  215. TMP = ABS(ZEIGSA(IZAMAX(M, ZEIGSA, 1))) ! The spectral radius of ZA
  216. ! Scale the matrix ZA to have unit spectral radius.
  217. CALL ZLASCL( 'G',0, 0, TMP, ONE, M, M, &
  218. ZA, LDA, INFO )
  219. CALL ZLASCL( 'G',0, 0, TMP, ONE, M, 1, &
  220. ZEIGSA, M, INFO )
  221. ANORM = ZLANGE( 'F', M, M, ZA, LDA, WDUMMY )
  222. IF ( K_TRAJ == 2 ) THEN
  223. ! generate data as two trajectories
  224. ! with two inital conditions
  225. CALL ZLARNV(2, ISEED, M, ZF(1,1) )
  226. DO i = 1, N/2
  227. CALL ZGEMV( 'N', M, M, ZONE, ZA, LDA, ZF(1,i), 1, &
  228. ZZERO, ZF(1,i+1), 1 )
  229. END DO
  230. ZX0(1:M,1:N/2) = ZF(1:M,1:N/2)
  231. ZY0(1:M,1:N/2) = ZF(1:M,2:N/2+1)
  232. CALL ZLARNV(2, ISEED, M, ZF(1,1) )
  233. DO i = 1, N-N/2
  234. CALL ZGEMV( 'N', M, M, ZONE, ZA, LDA, ZF(1,i), 1, &
  235. ZZERO, ZF(1,i+1), 1 )
  236. END DO
  237. ZX0(1:M,N/2+1:N) = ZF(1:M,1:N-N/2)
  238. ZY0(1:M,N/2+1:N) = ZF(1:M,2:N-N/2+1)
  239. ELSE
  240. CALL ZLARNV(2, ISEED, M, ZF(1,1) )
  241. DO i = 1, N
  242. CALL ZGEMV( 'N', M, M, ZONE, ZA, M, ZF(1,i), 1, &
  243. ZZERO, ZF(1,i+1), 1 )
  244. END DO
  245. ZF0(1:M,1:N+1) = ZF(1:M,1:N+1)
  246. ZX0(1:M,1:N) = ZF0(1:M,1:N)
  247. ZY0(1:M,1:N) = ZF0(1:M,2:N+1)
  248. END IF
  249. DEALLOCATE( ZEIGSA )
  250. !........................................................................
  251. DO iJOBZ = 1, 4
  252. SELECT CASE ( iJOBZ )
  253. CASE(1)
  254. JOBZ = 'V'
  255. RESIDS = 'R'
  256. CASE(2)
  257. JOBZ = 'V'
  258. RESIDS = 'N'
  259. CASE(3)
  260. JOBZ = 'F'
  261. RESIDS = 'N'
  262. CASE(4)
  263. JOBZ = 'N'
  264. RESIDS = 'N'
  265. END SELECT
  266. DO iJOBREF = 1, 3
  267. SELECT CASE ( iJOBREF )
  268. CASE(1)
  269. JOBREF = 'R'
  270. CASE(2)
  271. JOBREF = 'E'
  272. CASE(3)
  273. JOBREF = 'N'
  274. END SELECT
  275. DO iSCALE = 1, 4
  276. SELECT CASE ( iSCALE )
  277. CASE(1)
  278. SCALE = 'S'
  279. CASE(2)
  280. SCALE = 'C'
  281. CASE(3)
  282. SCALE = 'Y'
  283. CASE(4)
  284. SCALE = 'N'
  285. END SELECT
  286. DO iNRNK = -1, -2, -1
  287. NRNK = iNRNK
  288. NRNKsp = iNRNK
  289. DO iWHTSVD = 1, 3
  290. ! Check all four options to compute the POD basis
  291. ! via the SVD.
  292. WHTSVD = iWHTSVD
  293. WHTSVDsp = iWHTSVD
  294. DO LWMINOPT = 1, 2
  295. ! Workspace query for the minimal (1) and for the optimal
  296. ! (2) workspace lengths determined by workspace query.
  297. ! ZGEDMD is always tested and its results are also used for
  298. ! comparisons with ZGEDMDQ.
  299. ZX(1:M,1:N) = ZX0(1:M,1:N)
  300. ZY(1:M,1:N) = ZY0(1:M,1:N)
  301. CALL ZGEDMD( SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, &
  302. M, N, ZX, LDX, ZY, LDY, NRNK, TOL, &
  303. K, ZEIGS, ZZ, LDZ, RES, ZAU, LDAU, &
  304. ZW, LDW, ZS, LDS, ZDUMMY, -1, &
  305. WDUMMY, -1, IDUMMY, -1, INFO )
  306. IF ( (INFO .EQ. 2) .OR. ( INFO .EQ. 3 ) &
  307. .OR. ( INFO < 0 ) ) THEN
  308. WRITE(*,*) 'Call to ZGEDMD workspace query failed. &
  309. &Check the calling sequence and the code.'
  310. WRITE(*,*) 'The error code is ', INFO
  311. WRITE(*,*) 'The input parameters were ', &
  312. SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, &
  313. M, N, LDX, LDY, NRNK, TOL, LDZ, LDAU, LDW, LDS
  314. STOP
  315. END IF
  316. LZWORK = INT(ZDUMMY(LWMINOPT))
  317. LWORK = INT(WDUMMY(1))
  318. LIWORK = IDUMMY(1)
  319. ALLOCATE(ZWORK(LZWORK))
  320. ALLOCATE( WORK(LWORK))
  321. ALLOCATE(IWORK(LIWORK))
  322. CALL ZGEDMD( SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, &
  323. M, N, ZX, LDX, ZY, LDY, NRNK, TOL, &
  324. K, ZEIGS, ZZ, LDZ, RES, ZAU, LDAU, &
  325. ZW, LDW, ZS, LDS, ZWORK, LZWORK, &
  326. WORK, LWORK, IWORK, LIWORK, INFO )
  327. IF ( INFO /= 0 ) THEN
  328. WRITE(*,*) 'Call to ZGEDMD failed. &
  329. &Check the calling sequence and the code.'
  330. WRITE(*,*) 'The error code is ', INFO
  331. WRITE(*,*) 'The input parameters were ',&
  332. SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, &
  333. M, N, LDX, LDY, NRNK, TOL
  334. STOP
  335. END IF
  336. SINGVX(1:N) = WORK(1:N)
  337. !...... ZGEDMD check point
  338. IF ( LSAME(JOBZ,'V') ) THEN
  339. ! Check that Z = X*W, on return from ZGEDMD
  340. ! This checks that the returned eigenvectors in Z are
  341. ! the product of the SVD'POD basis returned in X
  342. ! and the eigenvectors of the rayleigh quotient
  343. ! returned in W
  344. CALL ZGEMM( 'N', 'N', M, K, K, ZONE, ZX, LDX, ZW, LDW, &
  345. ZZERO, ZZ1, LDZ )
  346. TMP = ZERO
  347. DO i = 1, K
  348. CALL ZAXPY( M, -ZONE, ZZ(1,i), 1, ZZ1(1,i), 1)
  349. TMP = MAX(TMP, DZNRM2( M, ZZ1(1,i), 1 ) )
  350. END DO
  351. TMP_ZXW = MAX(TMP_ZXW, TMP )
  352. IF ( TMP_ZXW <= 10*M*EPS ) THEN
  353. !WRITE(*,*) ' :) .... OK .........ZGEDMD PASSED.'
  354. ELSE
  355. NFAIL_Z_XV = NFAIL_Z_XV + 1
  356. WRITE(*,*) ':( .................ZGEDMD 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. !...... ZGEDMD 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 ZGEMM( 'N', 'N', M, K, M, ZONE, ZA, LDA, ZX, LDX, &
  373. ZZERO, ZZ1, LDZ )
  374. TMP = ZERO
  375. DO i = 1, K
  376. CALL ZAXPY( M, -ZONE, ZAU(1,i), 1, ZZ1(1,i), 1)
  377. TMP = MAX( TMP, DZNRM2( M, ZZ1(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. !WRITE(*,*) ':) .... OK .........ZGEDMD PASSED.'
  383. ELSE
  384. NFAIL_AU = NFAIL_AU + 1
  385. WRITE(*,*) ':( .................ZGEDMD FAILED!', &
  386. 'Check the code for implementation errors.'
  387. WRITE(*,*) 'The input parameters were ',&
  388. SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, &
  389. M, N, LDX, LDY, NRNK, TOL
  390. END IF
  391. ELSEIF ( LSAME(JOBREF,'E') ) THEN
  392. ! The unscaled vectors of the Exact DMD are computed.
  393. ! This option is included for the sake of completeness,
  394. ! for users who prefer the Exact DMD vectors. The
  395. ! returned vectors are in the real form, in the same way
  396. ! as the Ritz vectors. Here we just save the vectors
  397. ! and test them separately using a Matlab script.
  398. CALL ZGEMM( 'N', 'N', M, K, M, ZONE, ZA, LDA, ZAU, LDAU, ZZERO, ZY1, LDY )
  399. DO i=1, K
  400. ! have a real eigenvalue with real eigenvector
  401. CALL ZAXPY( M, -ZEIGS(i), ZAU(1,i), 1, ZY1(1,i), 1 )
  402. RESEX(i) = DZNRM2( M, ZY1(1,i), 1) / DZNRM2(M,ZAU(1,i),1)
  403. END DO
  404. END IF
  405. !...... ZGEDMD check point
  406. IF ( LSAME(RESIDS, 'R') ) THEN
  407. ! Compare the residuals returned by ZGEDMD with the
  408. ! explicitly computed residuals using the matrix A.
  409. ! Compute explicitly Y1 = A*Z
  410. CALL ZGEMM( 'N', 'N', M, K, M, ZONE, ZA, LDA, ZZ, LDZ, ZZERO, ZY1, LDY )
  411. ! ... and then A*Z(:,i) - LAMBDA(i)*Z(:,i), using the real forms
  412. ! of the invariant subspaces that correspond to complex conjugate
  413. ! pairs of eigencalues. (See the description of Z in ZGEDMD,)
  414. DO i=1, K
  415. ! have a real eigenvalue with real eigenvector
  416. CALL ZAXPY( M, -ZEIGS(i), ZZ(1,i), 1, ZY1(1,i), 1 )
  417. RES1(i) = DZNRM2( M, ZY1(1,i), 1)
  418. END DO
  419. TMP = ZERO
  420. DO i = 1, K
  421. TMP = MAX( TMP, ABS(RES(i) - RES1(i)) * &
  422. SINGVX(K)/(ANORM*SINGVX(1)) )
  423. END DO
  424. TMP_REZ = MAX( TMP_REZ, TMP )
  425. IF ( TMP <= TOL2 ) THEN
  426. !WRITE(*,*) ':) .... OK ..........ZGEDMD PASSED.'
  427. ELSE
  428. NFAIL_REZ = NFAIL_REZ + 1
  429. WRITE(*,*) ':( ..................ZGEDMD FAILED!', &
  430. 'Check the code for implementation errors.'
  431. WRITE(*,*) 'The input parameters were ',&
  432. SCALE, JOBZ, RESIDS, JOBREF, WHTSVD, &
  433. M, N, LDX, LDY, NRNK, TOL
  434. END IF
  435. IF ( LSAME(JOBREF,'E') ) THEN
  436. TMP = ZERO
  437. DO i = 1, K
  438. TMP = MAX( TMP, ABS(RES1(i) - RESEX(i))/(RES1(i)+RESEX(i)) )
  439. END DO
  440. TMP_EX = MAX(TMP_EX,TMP)
  441. END IF
  442. END IF
  443. DEALLOCATE(ZWORK)
  444. DEALLOCATE(WORK)
  445. DEALLOCATE(IWORK)
  446. IF ( TEST_QRDMD .AND. (K_TRAJ == 1) ) THEN
  447. ZF(1:M,1:N+1) = ZF0(1:M,1:N+1)
  448. CALL ZGEDMDQ( SCALE, JOBZ, RESIDS, WANTQ, WANTR, JOBREF, &
  449. WHTSVD, M, N+1, ZF, LDF, ZX, LDX, ZY, LDY, &
  450. NRNK, TOL, K, ZEIGS, ZZ, LDZ, RES, ZAU, &
  451. LDAU, ZW, LDW, ZS, LDS, ZDUMMY, -1, &
  452. WDUMMY, -1, IDUMMY, -1, INFO )
  453. LZWORK = INT(ZDUMMY(LWMINOPT))
  454. ALLOCATE( ZWORK(LZWORK) )
  455. LIWORK = IDUMMY(1)
  456. ALLOCATE(IWORK(LIWORK))
  457. LWORK = INT(WDUMMY(1))
  458. ALLOCATE(WORK(LWORK))
  459. CALL ZGEDMDQ( SCALE, JOBZ, RESIDS, WANTQ, WANTR, JOBREF, &
  460. WHTSVD, M, N+1, ZF, LDF, ZX, LDX, ZY, LDY, &
  461. NRNK, TOL, KQ, ZEIGS, ZZ, LDZ, RES, ZAU, &
  462. LDAU, ZW, LDW, ZS, LDS, ZWORK, LZWORK, &
  463. WORK, LWORK, IWORK, LIWORK, INFO )
  464. IF ( INFO /= 0 ) THEN
  465. WRITE(*,*) 'Call to ZGEDMDQ failed. &
  466. &Check the calling sequence and the code.'
  467. WRITE(*,*) 'The error code is ', INFO
  468. WRITE(*,*) 'The input parameters were ',&
  469. SCALE, JOBZ, RESIDS, WANTQ, WANTR, WHTSVD, &
  470. M, N, LDX, LDY, NRNK, TOL
  471. STOP
  472. END IF
  473. SINGVQX(1:N) = WORK(1:N)
  474. !..... ZGEDMDQ check point
  475. IF ( 1 == 0 ) THEN
  476. ! Comparison of ZGEDMD and ZGEDMDQ singular values disabled
  477. TMP = ZERO
  478. DO i = 1, MIN(K, KQ)
  479. TMP = MAX(TMP, ABS(SINGVX(i)-SINGVQX(i)) / &
  480. SINGVX(1) )
  481. END DO
  482. SVDIFF = MAX( SVDIFF, TMP )
  483. IF ( TMP > M*N*EPS ) THEN
  484. WRITE(*,*) 'FAILED! Something was wrong with the run.'
  485. NFAIL_SVDIFF = NFAIL_SVDIFF + 1
  486. DO j =1, 3
  487. write(*,*) j, SINGVX(j), SINGVQX(j)
  488. read(*,*)
  489. END DO
  490. END IF
  491. END IF
  492. !..... ZGEDMDQ check point
  493. IF ( LSAME(WANTQ,'Q') .AND. LSAME(WANTR,'R') ) THEN
  494. ! Check that the QR factors are computed and returned
  495. ! as requested. The residual ||F-Q*R||_F / ||F||_F
  496. ! is compared to M*N*EPS.
  497. ZF1(1:M,1:N+1) = ZF0(1:M,1:N+1)
  498. CALL ZGEMM( 'N', 'N', M, N+1, MIN(M,N+1), -ZONE, ZF, &
  499. LDF, ZY, LDY, ZONE, ZF1, LDF )
  500. TMP_FQR = ZLANGE( 'F', M, N+1, ZF1, LDF, WORK ) / &
  501. ZLANGE( 'F', M, N+1, ZF0, LDF, WORK )
  502. IF ( TMP_FQR > TOL2 ) THEN
  503. WRITE(*,*) 'FAILED! Something was wrong with the run.'
  504. NFAIL_F_QR = NFAIL_F_QR + 1
  505. ELSE
  506. !WRITE(*,*) '........ PASSED.'
  507. END IF
  508. END IF
  509. !..... ZGEDMDQ check point
  510. IF ( LSAME(RESIDS, 'R') ) THEN
  511. ! Compare the residuals returned by ZGEDMDQ with the
  512. ! explicitly computed residuals using the matrix A.
  513. ! Compute explicitly Y1 = A*Z
  514. CALL ZGEMM( 'N', 'N', M, KQ, M, ZONE, ZA, LDA, ZZ, LDZ, ZZERO, ZY1, LDY )
  515. ! ... and then A*Z(:,i) - LAMBDA(i)*Z(:,i), using the real forms
  516. ! of the invariant subspaces that correspond to complex conjugate
  517. ! pairs of eigencalues. (See the description of Z in ZGEDMDQ)
  518. DO i=1, KQ
  519. ! have a real eigenvalue with real eigenvector
  520. CALL ZAXPY( M, -ZEIGS(i), ZZ(1,i), 1, ZY1(1,i), 1 )
  521. ! Y(1:M,i) = Y(1:M,i) - REIG(i)*Z(1:M,i)
  522. RES1(i) = DZNRM2( M, ZY1(1,i), 1)
  523. END DO
  524. TMP = ZERO
  525. DO i = 1, KQ
  526. TMP = MAX( TMP, ABS(RES(i) - RES1(i)) * &
  527. SINGVQX(KQ)/(ANORM*SINGVQX(1)) )
  528. END DO
  529. TMP_REZQ = MAX( TMP_REZQ, TMP )
  530. IF ( TMP <= TOL2 ) THEN
  531. !WRITE(*,*) '.... OK ........ ZGEDMDQ PASSED.'
  532. ELSE
  533. NFAIL_REZQ = NFAIL_REZQ + 1
  534. WRITE(*,*) '................ ZGEDMDQ FAILED!', &
  535. 'Check the code for implementation errors.'
  536. STOP
  537. END IF
  538. END IF
  539. DEALLOCATE( ZWORK )
  540. DEALLOCATE( WORK )
  541. DEALLOCATE( IWORK )
  542. END IF ! ZGEDMDQ
  543. !.......................................................................................................
  544. END DO ! LWMINOPT
  545. !write(*,*) 'LWMINOPT loop completed'
  546. END DO ! iWHTSVD
  547. !write(*,*) 'WHTSVD loop completed'
  548. END DO ! iNRNK -2:-1
  549. !write(*,*) 'NRNK loop completed'
  550. END DO ! iSCALE 1:4
  551. !write(*,*) 'SCALE loop completed'
  552. END DO
  553. !write(*,*) 'JOBREF loop completed'
  554. END DO ! iJOBZ
  555. !write(*,*) 'JOBZ loop completed'
  556. END DO ! MODE -6:6
  557. !write(*,*) 'MODE loop completed'
  558. END DO ! 1 or 2 trajectories
  559. !write(*,*) 'trajectories loop completed'
  560. DEALLOCATE( ZA )
  561. DEALLOCATE( ZAC )
  562. DEALLOCATE( ZZ )
  563. DEALLOCATE( ZF )
  564. DEALLOCATE( ZF0 )
  565. DEALLOCATE( ZF1 )
  566. DEALLOCATE( ZX )
  567. DEALLOCATE( ZX0 )
  568. DEALLOCATE( ZY )
  569. DEALLOCATE( ZY0 )
  570. DEALLOCATE( ZY1 )
  571. DEALLOCATE( ZAU )
  572. DEALLOCATE( ZW )
  573. DEALLOCATE( ZS )
  574. DEALLOCATE( ZZ1 )
  575. DEALLOCATE( RES )
  576. DEALLOCATE( RES1 )
  577. DEALLOCATE( RESEX )
  578. DEALLOCATE( ZEIGS )
  579. DEALLOCATE( SINGVX )
  580. DEALLOCATE( SINGVQX )
  581. END DO ! LLOOP
  582. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  583. WRITE(*,*) ' Test summary for ZGEDMD :'
  584. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  585. WRITE(*,*)
  586. IF ( NFAIL_Z_XV == 0 ) THEN
  587. WRITE(*,*) '>>>> Z - U*V test PASSED.'
  588. ELSE
  589. WRITE(*,*) 'Z - U*V test FAILED ', NFAIL_Z_XV, ' time(s)'
  590. WRITE(*,*) 'Max error ||Z-U*V||_F was ', TMP_ZXW
  591. NFAIL_TOTAL = NFAIL_TOTAL + NFAIL_Z_XV
  592. END IF
  593. IF ( NFAIL_AU == 0 ) THEN
  594. WRITE(*,*) '>>>> A*U test PASSED. '
  595. ELSE
  596. WRITE(*,*) 'A*U test FAILED ', NFAIL_AU, ' time(s)'
  597. WRITE(*,*) 'Max A*U test adjusted error measure was ', TMP_AU
  598. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  599. NFAIL_TOTAL = NFAIL_TOTAL + NFAIL_AU
  600. END IF
  601. IF ( NFAIL_REZ == 0 ) THEN
  602. WRITE(*,*) '>>>> Rezidual computation test PASSED.'
  603. ELSE
  604. WRITE(*,*) 'Rezidual computation test FAILED ', NFAIL_REZ, 'time(s)'
  605. WRITE(*,*) 'Max residual computing test adjusted error measure was ', TMP_REZ
  606. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  607. NFAIL_TOTAL = NFAIL_TOTAL + NFAIL_REZ
  608. END IF
  609. IF ( NFAIL_TOTAL == 0 ) THEN
  610. WRITE(*,*) '>>>> ZGEDMD :: ALL TESTS PASSED.'
  611. ELSE
  612. WRITE(*,*) NFAIL_TOTAL, 'FAILURES!'
  613. WRITE(*,*) '>>>>>>>>>>>>>> ZGEDMD :: TESTS FAILED. CHECK THE IMPLEMENTATION.'
  614. END IF
  615. IF ( TEST_QRDMD ) THEN
  616. WRITE(*,*)
  617. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  618. WRITE(*,*) ' Test summary for ZGEDMDQ :'
  619. WRITE(*,*) '>>>>>>>>>>>>>>>>>>>>>>>>>>'
  620. WRITE(*,*)
  621. IF ( NFAIL_SVDIFF == 0 ) THEN
  622. WRITE(*,*) '>>>> ZGEDMD and ZGEDMDQ computed singular &
  623. &values test PASSED.'
  624. ELSE
  625. WRITE(*,*) 'ZGEDMD and ZGEDMDQ discrepancies in &
  626. &the singular values unacceptable ', &
  627. NFAIL_SVDIFF, ' times. Test FAILED.'
  628. WRITE(*,*) 'The maximal discrepancy in the singular values (relative to the norm) was ', SVDIFF
  629. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  630. NFAILQ_TOTAL = NFAILQ_TOTAL + NFAIL_SVDIFF
  631. END IF
  632. IF ( NFAIL_F_QR == 0 ) THEN
  633. WRITE(*,*) '>>>> F - Q*R test PASSED.'
  634. ELSE
  635. WRITE(*,*) 'F - Q*R test FAILED ', NFAIL_F_QR, ' time(s)'
  636. WRITE(*,*) 'The largest relative residual was ', TMP_FQR
  637. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  638. NFAILQ_TOTAL = NFAILQ_TOTAL + NFAIL_F_QR
  639. END IF
  640. IF ( NFAIL_REZQ == 0 ) THEN
  641. WRITE(*,*) '>>>> Rezidual computation test PASSED.'
  642. ELSE
  643. WRITE(*,*) 'Rezidual computation test FAILED ', NFAIL_REZQ, 'time(s)'
  644. WRITE(*,*) 'Max residual computing test adjusted error measure was ', TMP_REZQ
  645. WRITE(*,*) 'It should be up to O(M*N) times EPS, EPS = ', EPS
  646. NFAILQ_TOTAL = NFAILQ_TOTAL + NFAIL_REZQ
  647. END IF
  648. IF ( NFAILQ_TOTAL == 0 ) THEN
  649. WRITE(*,*) '>>>>>>> ZGEDMDQ :: ALL TESTS PASSED.'
  650. ELSE
  651. WRITE(*,*) NFAILQ_TOTAL, 'FAILURES!'
  652. WRITE(*,*) '>>>>>>> ZGEDMDQ :: TESTS FAILED. CHECK THE IMPLEMENTATION.'
  653. END IF
  654. END IF
  655. WRITE(*,*)
  656. WRITE(*,*) 'Test completed.'
  657. STOP
  658. END