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.

cchkdmd.f90 26 kB

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