493 SUBROUTINE zgedmd( JOBS, JOBZ, JOBR, JOBF, WHTSVD, &
494 M, N, X, LDX, Y, LDY, NRNK, TOL, &
495 K, EIGS, Z, LDZ, RES, B, LDB, &
496 W, LDW, S, LDS, ZWORK, LZWORK, &
497 RWORK, LRWORK, IWORK, LIWORK, INFO )
508 INTEGER,
PARAMETER :: WP = real64
512 CHARACTER,
INTENT(IN) :: JOBS, JOBZ, JOBR, JOBF
513 INTEGER,
INTENT(IN) :: WHTSVD, M, N, LDX, LDY, &
514 NRNK, LDZ, LDB, LDW, LDS, &
515 LIWORK, LRWORK, LZWORK
516 INTEGER,
INTENT(OUT) :: K, INFO
517 REAL(KIND=wp),
INTENT(IN) :: tol
521 COMPLEX(KIND=WP),
INTENT(INOUT) :: X(LDX,*), Y(LDY,*)
522 COMPLEX(KIND=WP),
INTENT(OUT) :: Z(LDZ,*), B(LDB,*), &
524 COMPLEX(KIND=WP),
INTENT(OUT) :: EIGS(*)
525 COMPLEX(KIND=WP),
INTENT(OUT) :: ZWORK(*)
526 REAL(KIND=wp),
INTENT(OUT) :: res(*)
527 REAL(KIND=wp),
INTENT(OUT) :: rwork(*)
528 INTEGER,
INTENT(OUT) :: IWORK(*)
532 REAL(KIND=wp),
PARAMETER :: one = 1.0_wp
533 REAL(KIND=wp),
PARAMETER :: zero = 0.0_wp
534 COMPLEX(KIND=WP),
PARAMETER :: ZONE = ( 1.0_wp, 0.0_wp )
535 COMPLEX(KIND=WP),
PARAMETER :: ZZERO = ( 0.0_wp, 0.0_wp )
539 REAL(KIND=wp) :: ofl, rootsc, scale, small, &
541 INTEGER :: i, j, IMINWR, INFO1, INFO2, &
542 LWRKEV, LWRSDD, LWRSVD, LWRSVJ, &
543 lwrsvq, mlwork, mwrkev, mwrsdd, &
544 mwrsvd, mwrsvj, mwrsvq, numrnk, &
546 LOGICAL :: BADXY, LQUERY, SCCOLX, SCCOLY, &
547 wntex, wntref, wntres, wntvec
548 CHARACTER :: JOBZL, T_OR_N
553 REAL(KIND=wp) :: rdummy(2)
560 LOGICAL DISNAN, LSAME
561 EXTERNAL disnan, lsame
571 INTRINSIC dble, int, max, sqrt
576 wntres = lsame(jobr,
'R')
577 sccolx = lsame(jobs,
'S') .OR. lsame(jobs,
'C')
578 sccoly = lsame(jobs,
'Y')
579 wntvec = lsame(jobz,
'V')
580 wntref = lsame(jobf,
'R')
581 wntex = lsame(jobf,
'E')
583 lquery = ( ( lzwork == -1 ) .OR. ( liwork == -1 ) &
584 .OR. ( lrwork == -1 ) )
586 IF ( .NOT. (sccolx .OR. sccoly .OR. &
587 lsame(jobs,
'N')) )
THEN
589 ELSE IF ( .NOT. (wntvec .OR. lsame(jobz,
'N') &
590 .OR. lsame(jobz,
'F')) )
THEN
592 ELSE IF ( .NOT. (wntres .OR. lsame(jobr,
'N')) .OR. &
593 ( wntres .AND. (.NOT.wntvec) ) )
THEN
595 ELSE IF ( .NOT. (wntref .OR. wntex .OR. &
596 lsame(jobf,
'N') ) )
THEN
598 ELSE IF ( .NOT.((whtsvd == 1) .OR. (whtsvd == 2) .OR. &
599 (whtsvd == 3) .OR. (whtsvd == 4) ))
THEN
601 ELSE IF ( m < 0 )
THEN
603 ELSE IF ( ( n < 0 ) .OR. ( n > m ) )
THEN
605 ELSE IF ( ldx < m )
THEN
607 ELSE IF ( ldy < m )
THEN
609 ELSE IF ( .NOT. (( nrnk == -2).OR.(nrnk == -1).OR. &
610 ((nrnk >= 1).AND.(nrnk <=n ))) )
THEN
612 ELSE IF ( ( tol < zero ) .OR. ( tol >= one ) )
THEN
614 ELSE IF ( ldz < m )
THEN
616 ELSE IF ( (wntref .OR. wntex ) .AND. ( ldb < m ) )
THEN
618 ELSE IF ( ldw < n )
THEN
620 ELSE IF ( lds < n )
THEN
624 IF ( info == 0 )
THEN
650 SELECT CASE ( whtsvd )
655 mwrsvd = max(1,2*min(m,n)+max(m,n))
656 mlwork = max(mlwork,mwrsvd)
657 mlrwrk = max(mlrwrk,n + 5*min(m,n))
659 CALL zgesvd(
'O',
'S', m, n, x, ldx, rwork, &
660 b, ldb, w, ldw, zwork, -1, rdummy, info1 )
661 lwrsvd = int( zwork(1) )
662 olwork = max(olwork,lwrsvd)
672 mwrsdd = 2*min(m,n)*min(m,n)+2*min(m,n)+max(m,n)
673 mlwork = max(mlwork,mwrsdd)
675 mlrwrk = max( mlrwrk, n + &
676 max( 5*min(m,n)*min(m,n)+7*min(m,n), &
677 5*min(m,n)*min(m,n)+5*min(m,n), &
678 2*max(m,n)*min(m,n)+ &
679 2*min(m,n)*min(m,n)+min(m,n) ) )
681 CALL zgesdd(
'O', m, n, x, ldx, rwork, b,ldb,&
682 w, ldw, zwork, -1, rdummy, iwork, info1 )
683 lwrsdd = max( mwrsdd,int( zwork(1) ))
685 olwork = max(olwork,lwrsdd)
688 CALL zgesvdq(
'H',
'P',
'N',
'R',
'R', m, n, &
689 x, ldx, rwork, z, ldz, w, ldw, numrnk, &
690 iwork, -1, zwork, -1, rdummy, -1, info1 )
692 mwrsvq = int(zwork(2))
693 mlwork = max(mlwork,mwrsvq)
694 mlrwrk = max(mlrwrk,n + int(rdummy(1)))
696 lwrsvq = int(zwork(1))
697 olwork = max(olwork,lwrsvq)
701 CALL zgejsv(
'F',
'U', jsvopt,
'R',
'N',
'P', m, &
702 n, x, ldx, rwork, z, ldz, w, ldw, &
703 zwork, -1, rdummy, -1, iwork, info1 )
705 mwrsvj = int(zwork(2))
706 mlwork = max(mlwork,mwrsvj)
707 mlrwrk = max(mlrwrk,n + max(7,int(rdummy(1))))
709 lwrsvj = int(zwork(1))
710 olwork = max(olwork,lwrsvj)
713 IF ( wntvec .OR. wntex .OR. lsame(jobz,
'F') )
THEN
719 mwrkev = max( 1, 2*n )
720 mlwork = max(mlwork,mwrkev)
721 mlrwrk = max(mlrwrk,n+2*n)
723 CALL zgeev(
'N', jobzl, n, s, lds, eigs, &
724 w, ldw, w, ldw, zwork, -1, rwork, info1 )
725 lwrkev = int(zwork(1))
726 olwork = max( olwork, lwrkev )
729 IF ( liwork < iminwr .AND. (.NOT.lquery) ) info = -30
730 IF ( lrwork < mlrwrk .AND. (.NOT.lquery) ) info = -28
731 IF ( lzwork < mlwork .AND. (.NOT.lquery) ) info = -26
736 CALL xerbla(
'ZGEDMD', -info )
738 ELSE IF ( lquery )
THEN
762 CALL zlassq( m, x(1,i), 1, scale, ssum )
763 IF ( disnan(scale) .OR. disnan(ssum) )
THEN
766 CALL xerbla(
'ZGEDMD',-info)
768 IF ( (scale /= zero) .AND. (ssum /= zero) )
THEN
770 IF ( scale .GE. (ofl / rootsc) )
THEN
781 CALL zlascl(
'G', 0, 0, scale, one/rootsc, &
782 m, 1, x(1,i), ldx, info2 )
783 rwork(i) = - scale * ( rootsc / dble(m) )
786 rwork(i) = scale * rootsc
787 CALL zlascl(
'G',0, 0, rwork(i), one, m, 1, &
801 CALL xerbla(
'ZGEDMD',-info)
806 IF ( rwork(i) > zero )
THEN
807 CALL zdscal( m, one/rwork(i), y(1,i), 1 )
809 ELSE IF ( rwork(i) < zero )
THEN
810 CALL zlascl(
'G', 0, 0, -rwork(i), &
811 one/dble(m), m, 1, y(1,i), ldy, info2 )
812 ELSE IF ( abs(y(izamax(m, y(1,i),1),i )) &
822 IF ( lsame(jobs,
'C')) &
823 CALL zdscal( m, zero, y(1,i), 1 )
835 CALL zlassq( m, y(1,i), 1, scale, ssum )
836 IF ( disnan(scale) .OR. disnan(ssum) )
THEN
839 CALL xerbla(
'ZGEDMD',-info)
841 IF ( scale /= zero .AND. (ssum /= zero) )
THEN
843 IF ( scale .GE. (ofl / rootsc) )
THEN
854 CALL zlascl(
'G', 0, 0, scale, one/rootsc, &
855 m, 1, y(1,i), ldy, info2 )
856 rwork(i) = - scale * ( rootsc / dble(m) )
859 rwork(i) = scale * rootsc
860 CALL zlascl(
'G',0, 0, rwork(i), one, m, 1, &
870 IF ( rwork(i) > zero )
THEN
871 CALL zdscal( m, one/rwork(i), x(1,i), 1 )
873 ELSE IF ( rwork(i) < zero )
THEN
874 CALL zlascl(
'G', 0, 0, -rwork(i), &
875 one/dble(m), m, 1, x(1,i), ldx, info2 )
876 ELSE IF ( abs(x(izamax(m, x(1,i),1),i )) &
894 SELECT CASE ( whtsvd )
896 CALL zgesvd(
'O',
'S', m, n, x, ldx, rwork, b, &
897 ldb, w, ldw, zwork, lzwork, rwork(n+1), info1 )
900 CALL zgesdd(
'O', m, n, x, ldx, rwork, b, ldb, w, &
901 ldw, zwork, lzwork, rwork(n+1), iwork, info1 )
904 CALL zgesvdq(
'H',
'P',
'N',
'R',
'R', m, n, &
905 x, ldx, rwork, z, ldz, w, ldw, &
906 numrnk, iwork, liwork, zwork, &
907 lzwork, rwork(n+1), lrwork-n, info1)
908 CALL zlacpy(
'A', m, numrnk, z, ldz, x, ldx )
911 CALL zgejsv(
'F',
'U', jsvopt,
'R',
'N',
'P', m, &
912 n, x, ldx, rwork, z, ldz, w, ldw, &
913 zwork, lzwork, rwork(n+1), lrwork-n, iwork, info1 )
914 CALL zlacpy(
'A', m, n, z, ldz, x, ldx )
918 IF ( xscl1 /= xscl2 )
THEN
925 CALL zlascl(
'G', 0, 0, xscl1, xscl2, m, n, y, ldy, info2
929 IF ( info1 > 0 )
THEN
936 IF ( rwork(1) == zero )
THEN
942 CALL xerbla(
'ZGEDMD',-info)
954 IF ( ( rwork(i) <= rwork(1)*tol ) .OR. &
955 ( rwork(i) <= small ) )
EXIT
961 IF ( ( rwork(i+1) <= rwork(i)*tol ) .OR. &
962 ( rwork(i) <= small ) )
EXIT
968 IF ( rwork(i) <= small )
EXIT
983 IF ( lsame(t_or_n,
'N') )
THEN
985 CALL zdscal( n, one/rwork(i), w(1,i), 1 )
997 rwork(n+i) = one/rwork(i)
1001 w(i,j) = cmplx(rwork(n+i),zero,kind=wp)*w(i,j)
1011 CALL zgemm(
'N', t_or_n, m, k, n, zone, y, ldy, w, &
1012 ldw, zzero, z, ldz )
1021 CALL zlacpy(
'A', m, k, z, ldz, b, ldb )
1024 CALL zgemm(
'C',
'N', k, k, m, zone, x, ldx, z, &
1025 ldz, zzero, s, lds )
1032 CALL zgemm(
'C',
'N', k, n, m, zone, x, ldx, y, ldy, &
1036 CALL zgemm(
'N', t_or_n, k, k, n, zone, z, ldz, w, &
1037 ldw, zzero, s, lds )
1043 IF ( wntres .OR. wntex )
THEN
1044 IF ( lsame(t_or_n,
'N') )
THEN
1045 CALL zlacpy(
'A', n, k, w, ldw, z, ldz )
1047 CALL zlacpy(
'A', k, n, w, ldw, z, ldz )
1055 CALL zgeev(
'N', jobzl, k, s, lds, eigs, w, ldw, &
1056 w, ldw, zwork, lzwork, rwork(n+1), info1 )
1061 IF ( info1 > 0 )
THEN
1071 IF ( wntvec .OR. wntex )
THEN
1077 CALL zgemm(
'N',
'N', m, k, k, zone, z, ldz, w, &
1078 ldw, zzero, y, ldy )
1084 CALL zgemm( t_or_n,
'N', n, k, k, zone, z, ldz, &
1085 w, ldw, zzero, s, lds )
1089 CALL zgemm(
'N',
'N', m, k, n, zone, y, ldy, s, &
1090 lds, zzero, z, ldz )
1093 CALL zlacpy(
'A', m, k, z, ldz, y, ldy )
1094 IF ( wntex )
CALL zlacpy(
'A', m, k, z, ldz, b, ldb )
1096 ELSE IF ( wntex )
THEN
1099 CALL zgemm( t_or_n,
'N', n, k, k, zone, z, ldz, &
1100 w, ldw, zzero, s, lds )
1104 CALL zgemm(
'N',
'N', m, k, n, zone, y, ldy, s, &
1105 lds, zzero, b, ldb )
1116 IF ( wntvec )
CALL zgemm(
'N',
'N', m, k, k, zone, x, ldx, w, ldw,
1122 CALL zaxpy( m, -eigs(i), z(1,i), 1, y(1,i), 1 )
1124 res(i) =
dznrm2( m, y(1,i), 1 )
1129 IF ( whtsvd == 4 )
THEN
1135 IF ( .NOT. badxy )
THEN
subroutine xerbla(srname, info)
subroutine zaxpy(n, za, zx, incx, zy, incy)
ZAXPY
subroutine zgedmd(jobs, jobz, jobr, jobf, whtsvd, m, n, x, ldx, y, ldy, nrnk, tol, k, eigs, z, ldz, res, b, ldb, w, ldw, s, lds, zwork, lzwork, rwork, lrwork, iwork, liwork, info)
ZGEDMD computes the Dynamic Mode Decomposition (DMD) for a pair of data snapshot matrices.
subroutine zgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
ZGEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices
subroutine zgejsv(joba, jobu, jobv, jobr, jobt, jobp, m, n, a, lda, sva, u, ldu, v, ldv, cwork, lwork, rwork, lrwork, iwork, info)
ZGEJSV
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
subroutine zgesdd(jobz, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, rwork, iwork, info)
ZGESDD
subroutine zgesvd(jobu, jobvt, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, rwork, info)
ZGESVD computes the singular value decomposition (SVD) for GE matrices
subroutine zgesvdq(joba, jobp, jobr, jobu, jobv, m, n, a, lda, s, u, ldu, v, ldv, numrank, iwork, liwork, cwork, lcwork, rwork, lrwork, info)
ZGESVDQ computes the singular value decomposition (SVD) with a QR-Preconditioned QR SVD Method for GE...
integer function izamax(n, zx, incx)
IZAMAX
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
double precision function dlamch(cmach)
DLAMCH
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
subroutine zlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
ZLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
subroutine zlassq(n, x, incx, scale, sumsq)
ZLASSQ updates a sum of squares represented in scaled form.
real(wp) function dznrm2(n, x, incx)
DZNRM2
subroutine zdscal(n, da, zx, incx)
ZDSCAL