496 SUBROUTINE cgedmd( JOBS, JOBZ, JOBR, JOBF, WHTSVD, &
497 M, N, X, LDX, Y, LDY, NRNK, TOL, &
498 K, EIGS, Z, LDZ, RES, B, LDB, &
499 W, LDW, S, LDS, ZWORK, LZWORK, &
500 RWORK, LRWORK, IWORK, LIWORK, INFO )
511 INTEGER,
PARAMETER :: WP = real32
515 CHARACTER,
INTENT(IN) :: JOBS, JOBZ, JOBR, JOBF
516 INTEGER,
INTENT(IN) :: WHTSVD, M, N, LDX, LDY, &
517 NRNK, LDZ, LDB, LDW, LDS, &
518 LIWORK, LRWORK, LZWORK
519 INTEGER,
INTENT(OUT) :: K, INFO
520 REAL(KIND=wp),
INTENT(IN) :: tol
524 COMPLEX(KIND=WP),
INTENT(INOUT) :: X(LDX,*), Y(LDY,*)
525 COMPLEX(KIND=WP),
INTENT(OUT) :: Z(LDZ,*), B(LDB,*), &
527 COMPLEX(KIND=WP),
INTENT(OUT) :: EIGS(*)
528 COMPLEX(KIND=WP),
INTENT(OUT) :: ZWORK(*)
529 REAL(KIND=wp),
INTENT(OUT) :: res(*)
530 REAL(KIND=wp),
INTENT(OUT) :: rwork(*)
531 INTEGER,
INTENT(OUT) :: IWORK(*)
535 REAL(KIND=wp),
PARAMETER :: one = 1.0_wp
536 REAL(KIND=wp),
PARAMETER :: zero = 0.0_wp
537 COMPLEX(KIND=WP),
PARAMETER :: ZONE = ( 1.0_wp, 0.0_wp )
538 COMPLEX(KIND=WP),
PARAMETER :: ZZERO = ( 0.0_wp, 0.0_wp )
542 REAL(KIND=wp) :: ofl, rootsc, scale, small, &
544 INTEGER :: i, j, IMINWR, INFO1, INFO2, &
545 LWRKEV, LWRSDD, LWRSVD, LWRSVJ, &
546 lwrsvq, mlwork, mwrkev, mwrsdd, &
547 mwrsvd, mwrsvj, mwrsvq, numrnk, &
549 LOGICAL :: BADXY, LQUERY, SCCOLX, SCCOLY, &
550 wntex, wntref, wntres, wntvec
551 CHARACTER :: JOBZL, T_OR_N
556 REAL(KIND=wp) :: rdummy(2)
563 LOGICAL SISNAN, LSAME
564 EXTERNAL sisnan, lsame
574 INTRINSIC float, int, max, sqrt
579 wntres = lsame(jobr,
'R')
580 sccolx = lsame(jobs,
'S') .OR. lsame(jobs,
'C')
581 sccoly = lsame(jobs,
'Y')
582 wntvec = lsame(jobz,
'V')
583 wntref = lsame(jobf,
'R')
584 wntex = lsame(jobf,
'E')
586 lquery = ( ( lzwork == -1 ) .OR. ( liwork == -1 ) &
587 .OR. ( lrwork == -1 ) )
589 IF ( .NOT. (sccolx .OR. sccoly .OR. &
590 lsame(jobs,
'N')) )
THEN
592 ELSE IF ( .NOT. (wntvec .OR. lsame(jobz,
'N') &
593 .OR. lsame(jobz,
'F')) )
THEN
595 ELSE IF ( .NOT. (wntres .OR. lsame(jobr,
'N')) .OR. &
596 ( wntres .AND. (.NOT.wntvec) ) )
THEN
598 ELSE IF ( .NOT. (wntref .OR. wntex .OR. &
599 lsame(jobf,
'N') ) )
THEN
601 ELSE IF ( .NOT.((whtsvd == 1) .OR. (whtsvd == 2) .OR. &
602 (whtsvd == 3) .OR. (whtsvd == 4) ))
THEN
604 ELSE IF ( m < 0 )
THEN
606 ELSE IF ( ( n < 0 ) .OR. ( n > m ) )
THEN
608 ELSE IF ( ldx < m )
THEN
610 ELSE IF ( ldy < m )
THEN
612 ELSE IF ( .NOT. (( nrnk == -2).OR.(nrnk == -1).OR. &
613 ((nrnk >= 1).AND.(nrnk <=n ))) )
THEN
615 ELSE IF ( ( tol < zero ) .OR. ( tol >= one ) )
THEN
617 ELSE IF ( ldz < m )
THEN
619 ELSE IF ( (wntref .OR. wntex ) .AND. ( ldb < m ) )
THEN
621 ELSE IF ( ldw < n )
THEN
623 ELSE IF ( lds < n )
THEN
627 IF ( info == 0 )
THEN
653 SELECT CASE ( whtsvd )
658 mwrsvd = max(1,2*min(m,n)+max(m,n))
659 mlwork = max(mlwork,mwrsvd)
660 mlrwrk = max(mlrwrk,n + 5*min(m,n))
662 CALL cgesvd(
'O',
'S', m, n, x, ldx, rwork, &
663 b, ldb, w, ldw, zwork, -1, rdummy, info1 )
664 lwrsvd = int( zwork(1) )
665 olwork = max(olwork,lwrsvd)
675 mwrsdd = 2*min(m,n)*min(m,n)+2*min(m,n)+max(m,n)
676 mlwork = max(mlwork,mwrsdd)
678 mlrwrk = max( mlrwrk, n + &
679 max( 5*min(m,n)*min(m,n)+7*min(m,n), &
680 5*min(m,n)*min(m,n)+5*min(m,n), &
681 2*max(m,n)*min(m,n)+ &
682 2*min(m,n)*min(m,n)+min(m,n) ) )
684 CALL cgesdd(
'O', m, n, x, ldx, rwork, b, &
685 ldb, w, ldw, zwork, -1, rdummy, iwork, info1 )
686 lwrsdd = max(mwrsdd,int( zwork(1) ))
687 olwork = max(olwork,lwrsdd)
690 CALL cgesvdq(
'H',
'P',
'N',
'R',
'R', m, n, &
691 x, ldx, rwork, z, ldz, w, ldw, numrnk, &
692 iwork, -1, zwork, -1, rdummy, -1, info1 )
694 mwrsvq = int(zwork(2))
695 mlwork = max(mlwork,mwrsvq)
696 mlrwrk = max(mlrwrk,n + int(rdummy(1)))
698 lwrsvq = int(zwork(1))
699 olwork = max(olwork,lwrsvq)
703 CALL cgejsv(
'F',
'U', jsvopt,
'N',
'N',
'P', m, &
704 n, x, ldx, rwork, z, ldz, w, ldw, &
705 zwork, -1, rdummy, -1, iwork, info1 )
707 mwrsvj = int(zwork(2))
708 mlwork = max(mlwork,mwrsvj)
709 mlrwrk = max(mlrwrk,n + max(7,int(rdummy(1))))
711 lwrsvj = int(zwork(1))
712 olwork = max(olwork,lwrsvj)
715 IF ( wntvec .OR. wntex .OR. lsame(jobz,
'F') )
THEN
721 mwrkev = max( 1, 2*n )
722 mlwork = max(mlwork,mwrkev)
723 mlrwrk = max(mlrwrk,n+2*n)
725 CALL cgeev(
'N', jobzl, n, s, lds, eigs, &
726 w, ldw, w, ldw, zwork, -1, rwork, info1 )
727 lwrkev = int(zwork(1))
728 olwork = max( olwork, lwrkev )
729 olwork = max( 2, olwork )
732 IF ( liwork < iminwr .AND. (.NOT.lquery) ) info = -30
733 IF ( lrwork < mlrwrk .AND. (.NOT.lquery) ) info = -28
734 IF ( lzwork < mlwork .AND. (.NOT.lquery) ) info = -26
739 CALL xerbla(
'CGEDMD', -info )
741 ELSE IF ( lquery )
THEN
765 CALL classq( m, x(1,i), 1, scale, ssum )
766 IF ( sisnan(scale) .OR. sisnan(ssum) )
THEN
769 CALL xerbla(
'CGEDMD',-info)
771 IF ( (scale /= zero) .AND. (ssum /= zero) )
THEN
773 IF ( scale .GE. (ofl / rootsc) )
THEN
784 CALL clascl(
'G', 0, 0, scale, one/rootsc, &
785 m, 1, x(1,i), ldx, info2 )
786 rwork(i) = - scale * ( rootsc / float(m) )
789 rwork(i) = scale * rootsc
790 CALL clascl(
'G',0, 0, rwork(i), one, m, 1, &
804 CALL xerbla(
'CGEDMD',-info)
809 IF ( rwork(i) > zero )
THEN
810 CALL csscal( m, one/rwork(i), y(1,i), 1 )
812 ELSE IF ( rwork(i) < zero )
THEN
813 CALL clascl(
'G', 0, 0, -rwork(i), &
814 one/float(m), m, 1, y(1,i), ldy, info2 )
815 ELSE IF ( abs(y(icamax(m, y(1,i),1),i )) &
825 IF ( lsame(jobs,
'C')) &
826 CALL csscal( m, zero, y(1,i), 1 )
838 CALL classq( m, y(1,i), 1, scale, ssum )
839 IF ( sisnan(scale) .OR. sisnan(ssum) )
THEN
842 CALL xerbla(
'CGEDMD',-info)
844 IF ( scale /= zero .AND. (ssum /= zero) )
THEN
846 IF ( scale .GE. (ofl / rootsc) )
THEN
857 CALL clascl(
'G', 0, 0, scale, one/rootsc, &
858 m, 1, y(1,i), ldy, info2 )
859 rwork(i) = - scale * ( rootsc / float(m) )
862 rwork(i) = scale * rootsc
863 CALL clascl(
'G',0, 0, rwork(i), one, m, 1, &
873 IF ( rwork(i) > zero )
THEN
874 CALL csscal( m, one/rwork(i), x(1,i), 1 )
876 ELSE IF ( rwork(i) < zero )
THEN
877 CALL clascl(
'G', 0, 0, -rwork(i), &
878 one/float(m), m, 1, x(1,i), ldx, info2 )
879 ELSE IF ( abs(x(icamax(m, x(1,i),1),i )) &
897 SELECT CASE ( whtsvd )
899 CALL cgesvd(
'O',
'S', m, n, x, ldx, rwork, b, &
900 ldb, w, ldw, zwork, lzwork, rwork(n+1), info1 )
903 CALL cgesdd(
'O', m, n, x, ldx, rwork, b, ldb, w, &
904 ldw, zwork, lzwork, rwork(n+1), iwork, info1 )
907 CALL cgesvdq(
'H',
'P',
'N',
'R',
'R', m, n, &
908 x, ldx, rwork, z, ldz, w, ldw, &
909 numrnk, iwork, liwork, zwork, &
910 lzwork, rwork(n+1), lrwork-n, info1)
911 CALL clacpy(
'A', m, numrnk, z, ldz, x, ldx )
914 CALL cgejsv(
'F',
'U', jsvopt,
'N',
'N',
'P', m, &
915 n, x, ldx, rwork, z, ldz, w, ldw, &
916 zwork, lzwork, rwork(n+1), lrwork-n, iwork, info1 )
917 CALL clacpy(
'A', m, n, z, ldz, x, ldx )
921 IF ( xscl1 /= xscl2 )
THEN
928 CALL clascl(
'G', 0, 0, xscl1, xscl2, m, n, y, ldy, info2
932 IF ( info1 > 0 )
THEN
939 IF ( rwork(1) == zero )
THEN
945 CALL xerbla(
'CGEDMD',-info)
957 IF ( ( rwork(i) <= rwork(1)*tol ) .OR. &
958 ( rwork(i) <= small ) )
EXIT
964 IF ( ( rwork(i+1) <= rwork(i)*tol ) .OR. &
965 ( rwork(i) <= small ) )
EXIT
971 IF ( rwork(i) <= small )
EXIT
986 IF ( lsame(t_or_n,
'N') )
THEN
988 CALL csscal( n, one/rwork(i), w(1,i), 1 )
1000 rwork(n+i) = one/rwork(i)
1004 w(i,j) = cmplx(rwork(n+i),zero,kind=wp)*w(i,j)
1014 CALL cgemm(
'N', t_or_n, m, k, n, zone, y, ldy, w, &
1015 ldw, zzero, z, ldz )
1024 CALL clacpy(
'A', m, k, z, ldz, b, ldb )
1027 CALL cgemm(
'C',
'N', k, k, m, zone, x, ldx, z, &
1028 ldz, zzero, s, lds )
1035 CALL cgemm(
'C',
'N', k, n, m, zone, x, ldx, y, ldy, &
1039 CALL cgemm(
'N', t_or_n, k, k, n, zone, z, ldz, w, &
1040 ldw, zzero, s, lds )
1046 IF ( wntres .OR. wntex )
THEN
1047 IF ( lsame(t_or_n,
'N') )
THEN
1048 CALL clacpy(
'A', n, k, w, ldw, z, ldz )
1050 CALL clacpy(
'A', k, n, w, ldw, z, ldz )
1058 CALL cgeev(
'N', jobzl, k, s, lds, eigs, w, &
1059 ldw, w, ldw, zwork, lzwork, rwork(n+1), info1 )
1064 IF ( info1 > 0 )
THEN
1074 IF ( wntvec .OR. wntex )
THEN
1080 CALL cgemm(
'N',
'N', m, k, k, zone, z, ldz, w, &
1081 ldw, zzero, y, ldy )
1087 CALL cgemm( t_or_n,
'N', n, k, k, zone, z, ldz, &
1088 w, ldw, zzero, s, lds)
1092 CALL cgemm(
'N',
'N', m, k, n, zone, y, ldy, s, &
1096 CALL clacpy(
'A', m, k, z, ldz, y, ldy )
1097 IF ( wntex )
CALL clacpy(
'A', m, k, z, ldz, b, ldb )
1099 ELSE IF ( wntex )
THEN
1102 CALL cgemm( t_or_n,
'N', n, k, k, zone, z, ldz, &
1103 w, ldw, zzero, s, lds)
1107 CALL cgemm(
'N',
'N', m, k, n, zone, y, ldy, s, &
1119 IF ( wntvec )
CALL cgemm(
'N',
'N', m, k, k, zone, x, ldx, w, ldw
1125 CALL caxpy( m, -eigs(i), z(1,i), 1, y(1,i), 1 )
1127 res(i) =
scnrm2( m, y(1,i), 1)
1132 IF ( whtsvd == 4 )
THEN
1138 IF ( .NOT. badxy )
THEN
subroutine xerbla(srname, info)
subroutine caxpy(n, ca, cx, incx, cy, incy)
CAXPY
subroutine cgedmd(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)
CGEDMD computes the Dynamic Mode Decomposition (DMD) for a pair of data snapshot matrices.
subroutine cgeev(jobvl, jobvr, n, a, lda, w, vl, ldvl, vr, ldvr, work, lwork, rwork, info)
CGEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices
subroutine cgejsv(joba, jobu, jobv, jobr, jobt, jobp, m, n, a, lda, sva, u, ldu, v, ldv, cwork, lwork, rwork, lrwork, iwork, info)
CGEJSV
subroutine cgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CGEMM
subroutine cgesdd(jobz, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, rwork, iwork, info)
CGESDD
subroutine cgesvd(jobu, jobvt, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, rwork, info)
CGESVD computes the singular value decomposition (SVD) for GE matrices
subroutine cgesvdq(joba, jobp, jobr, jobu, jobv, m, n, a, lda, s, u, ldu, v, ldv, numrank, iwork, liwork, cwork, lcwork, rwork, lrwork, info)
CGESVDQ computes the singular value decomposition (SVD) with a QR-Preconditioned QR SVD Method for GE...
integer function icamax(n, cx, incx)
ICAMAX
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
real function slamch(cmach)
SLAMCH
real function clange(norm, m, n, a, lda, work)
CLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
subroutine clascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
CLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
subroutine classq(n, x, incx, scale, sumsq)
CLASSQ updates a sum of squares represented in scaled form.
real(wp) function scnrm2(n, x, incx)
SCNRM2
subroutine csscal(n, sa, cx, incx)
CSSCAL