LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sgedmd.f90
Go to the documentation of this file.
1
2!
3! =========== DOCUMENTATION ===========
4!
5! Definition:
6! ===========
7!
8! SUBROUTINE SGEDMD( JOBS, JOBZ, JOBR, JOBF, WHTSVD, &
9! M, N, X, LDX, Y, LDY, NRNK, TOL, &
10! K, REIG, IMEIG, Z, LDZ, RES, &
11! B, LDB, W, LDW, S, LDS, &
12! WORK, LWORK, IWORK, LIWORK, INFO )
13!.....
14! USE iso_fortran_env
15! IMPLICIT NONE
16! INTEGER, PARAMETER :: WP = real32
17!.....
18! Scalar arguments
19! CHARACTER, INTENT(IN) :: JOBS, JOBZ, JOBR, JOBF
20! INTEGER, INTENT(IN) :: WHTSVD, M, N, LDX, LDY, &
21! NRNK, LDZ, LDB, LDW, LDS, &
22! LWORK, LIWORK
23! INTEGER, INTENT(OUT) :: K, INFO
24! REAL(KIND=WP), INTENT(IN) :: TOL
25! Array arguments
26! REAL(KIND=WP), INTENT(INOUT) :: X(LDX,*), Y(LDY,*)
27! REAL(KIND=WP), INTENT(OUT) :: Z(LDZ,*), B(LDB,*), &
28! W(LDW,*), S(LDS,*)
29! REAL(KIND=WP), INTENT(OUT) :: REIG(*), IMEIG(*), &
30! RES(*)
31! REAL(KIND=WP), INTENT(OUT) :: WORK(*)
32! INTEGER, INTENT(OUT) :: IWORK(*)
33!
34!............................................................
36! =============
51!............................................................
53! ================
69!......................................................................
71! ================================
91!......................................................................
93! ==============================
99!============================================================
100! Arguments
101! =========
102!
121!.....
138!.....
149!.....
163!.....
187!.....
193!.....
200!.....
213!.....
219!.....
231!.....
237!.....
262!.....
269!.....
279!.....
288!.....
307!.....
335!.....
341!.....
361!.....
374!.....
380!.....
393!.....
399!.....
408!.....
414!.....
432!.....
475!.....
485!.....
499!.....
520!
521! Authors:
522! ========
523!
525!
527!
528!.............................................................
529!.............................................................
530 SUBROUTINE sgedmd( JOBS, JOBZ, JOBR, JOBF, WHTSVD, &
531 M, N, X, LDX, Y, LDY, NRNK, TOL, &
532 K, REIG, IMEIG, Z, LDZ, RES, &
533 B, LDB, W, LDW, S, LDS, &
534 WORK, LWORK, IWORK, LIWORK, INFO )
535!
536! -- LAPACK driver routine --
537!
538! -- LAPACK is a software package provided by University of --
539! -- Tennessee, University of California Berkeley, University of --
540! -- Colorado Denver and NAG Ltd.. --
541!
542!.....
543 USE iso_fortran_env
544 IMPLICIT NONE
545 INTEGER, PARAMETER :: WP = real32
546!
547! Scalar arguments
548! ~~~~~~~~~~~~~~~~
549 CHARACTER, INTENT(IN) :: JOBS, JOBZ, JOBR, JOBF
550 INTEGER, INTENT(IN) :: WHTSVD, M, N, LDX, LDY, &
551 NRNK, LDZ, LDB, LDW, LDS, &
552 LWORK, LIWORK
553 INTEGER, INTENT(OUT) :: K, INFO
554 REAL(KIND=wp), INTENT(IN) :: tol
555!
556! Array arguments
557! ~~~~~~~~~~~~~~~
558 REAL(KIND=wp), INTENT(INOUT) :: x(ldx,*), y(ldy,*)
559 REAL(KIND=wp), INTENT(OUT) :: z(ldz,*), b(ldb,*), &
560 w(ldw,*), s(lds,*)
561 REAL(KIND=wp), INTENT(OUT) :: reig(*), imeig(*), &
562 res(*)
563 REAL(KIND=wp), INTENT(OUT) :: work(*)
564 INTEGER, INTENT(OUT) :: IWORK(*)
565!
566! Parameters
567! ~~~~~~~~~~
568 REAL(KIND=wp), PARAMETER :: one = 1.0_wp
569 REAL(KIND=wp), PARAMETER :: zero = 0.0_wp
570!
571! Local scalars
572! ~~~~~~~~~~~~~
573 REAL(KIND=wp) :: ofl, rootsc, scale, small, &
574 ssum, xscl1, xscl2
575 INTEGER :: i, j, IMINWR, INFO1, INFO2, &
576 LWRKEV, LWRSDD, LWRSVD, &
577 lwrsvq, mlwork, mwrkev, mwrsdd, &
578 mwrsvd, mwrsvj, mwrsvq, numrnk, &
579 olwork
580 LOGICAL :: BADXY, LQUERY, SCCOLX, SCCOLY, &
581 wntex, wntref, wntres, wntvec
582 CHARACTER :: JOBZL, T_OR_N
583 CHARACTER :: JSVOPT
584!
585! Local arrays
586! ~~~~~~~~~~~~
587 REAL(KIND=wp) :: ab(2,2), rdummy(2), rdummy2(2)
588!
589! External functions (BLAS and LAPACK)
590! ~~~~~~~~~~~~~~~~~
591 REAL(KIND=wp) slange, slamch, snrm2
592 EXTERNAL slange, slamch, snrm2, isamax
593 INTEGER ISAMAX
594 LOGICAL SISNAN, LSAME
595 EXTERNAL sisnan, lsame
596!
597! External subroutines (BLAS and LAPACK)
598! ~~~~~~~~~~~~~~~~~~~~
599 EXTERNAL saxpy, sgemm, sscal
600 EXTERNAL sgeev, sgejsv, sgesdd, sgesvd, sgesvdq, &
602!
603! Intrinsic functions
604! ~~~~~~~~~~~~~~~~~~~
605 INTRINSIC int, float, max, sqrt
606!............................................................
607!
608! Test the input arguments
609!
610 wntres = lsame(jobr,'R')
611 sccolx = lsame(jobs,'S') .OR. lsame(jobs,'C')
612 sccoly = lsame(jobs,'Y')
613 wntvec = lsame(jobz,'V')
614 wntref = lsame(jobf,'R')
615 wntex = lsame(jobf,'E')
616 info = 0
617 lquery = ( ( lwork == -1 ) .OR. ( liwork == -1 ) )
618!
619 IF ( .NOT. (sccolx .OR. sccoly .OR. &
620 lsame(jobs,'N')) ) THEN
621 info = -1
622 ELSE IF ( .NOT. (wntvec .OR. lsame(jobz,'N') &
623 .OR. lsame(jobz,'F')) ) THEN
624 info = -2
625 ELSE IF ( .NOT. (wntres .OR. lsame(jobr,'N')) .OR. &
626 ( wntres .AND. (.NOT.wntvec) ) ) THEN
627 info = -3
628 ELSE IF ( .NOT. (wntref .OR. wntex .OR. &
629 lsame(jobf,'N') ) ) THEN
630 info = -4
631 ELSE IF ( .NOT.((whtsvd == 1) .OR. (whtsvd == 2) .OR. &
632 (whtsvd == 3) .OR. (whtsvd == 4) )) THEN
633 info = -5
634 ELSE IF ( m < 0 ) THEN
635 info = -6
636 ELSE IF ( ( n < 0 ) .OR. ( n > m ) ) THEN
637 info = -7
638 ELSE IF ( ldx < m ) THEN
639 info = -9
640 ELSE IF ( ldy < m ) THEN
641 info = -11
642 ELSE IF ( .NOT. (( nrnk == -2).OR.(nrnk == -1).OR. &
643 ((nrnk >= 1).AND.(nrnk <=n ))) ) THEN
644 info = -12
645 ELSE IF ( ( tol < zero ) .OR. ( tol >= one ) ) THEN
646 info = -13
647 ELSE IF ( ldz < m ) THEN
648 info = -18
649 ELSE IF ( (wntref .OR. wntex ) .AND. ( ldb < m ) ) THEN
650 info = -21
651 ELSE IF ( ldw < n ) THEN
652 info = -23
653 ELSE IF ( lds < n ) THEN
654 info = -25
655 END IF
656!
657 IF ( info == 0 ) THEN
658 ! Compute the minimal and the optimal workspace
659 ! requirements. Simulate running the code and
660 ! determine minimal and optimal sizes of the
661 ! workspace at any moment of the run.
662 IF ( n == 0 ) THEN
663 ! Quick return. All output except K is void.
664 ! INFO=1 signals the void input.
665 ! In case of a workspace query, the default
666 ! minimal workspace lengths are returned.
667 IF ( lquery ) THEN
668 iwork(1) = 1
669 work(1) = 2
670 work(2) = 2
671 ELSE
672 k = 0
673 END IF
674 info = 1
675 RETURN
676 END IF
677 mlwork = max(2,n)
678 olwork = max(2,n)
679 iminwr = 1
680 SELECT CASE ( whtsvd )
681 CASE (1)
682 ! The following is specified as the minimal
683 ! length of WORK in the definition of SGESVD:
684 ! MWRSVD = MAX(1,3*MIN(M,N)+MAX(M,N),5*MIN(M,N))
685 mwrsvd = max(1,3*min(m,n)+max(m,n),5*min(m,n))
686 mlwork = max(mlwork,n + mwrsvd)
687 IF ( lquery ) THEN
688 CALL sgesvd( 'O', 'S', m, n, x, ldx, work, &
689 b, ldb, w, ldw, rdummy, -1, info1 )
690 lwrsvd = max( mwrsvd, int( rdummy(1) ) )
691 olwork = max(olwork,n + lwrsvd)
692 END IF
693 CASE (2)
694 ! The following is specified as the minimal
695 ! length of WORK in the definition of SGESDD:
696 ! MWRSDD = 3*MIN(M,N)*MIN(M,N) +
697 ! MAX( MAX(M,N),5*MIN(M,N)*MIN(M,N)+4*MIN(M,N) )
698 ! IMINWR = 8*MIN(M,N)
699 mwrsdd = 3*min(m,n)*min(m,n) + &
700 max( max(m,n),5*min(m,n)*min(m,n)+4*min(m,n) )
701 mlwork = max(mlwork,n + mwrsdd)
702 iminwr = 8*min(m,n)
703 IF ( lquery ) THEN
704 CALL sgesdd( 'O', m, n, x, ldx, work, b, &
705 ldb, w, ldw, rdummy, -1, iwork, info1 )
706 lwrsdd = max( mwrsdd, int( rdummy(1) ) )
707 olwork = max(olwork,n + lwrsdd)
708 END IF
709 CASE (3)
710 !LWQP3 = 3*N+1
711 !LWORQ = MAX(N, 1)
712 !MWRSVD = MAX(1,3*MIN(M,N)+MAX(M,N),5*MIN(M,N))
713 !MWRSVQ = N + MAX( LWQP3, MWRSVD, LWORQ )+ MAX(M,2)
714 !MLWORK = N + MWRSVQ
715 !IMINWR = M+N-1
716 CALL sgesvdq( 'H', 'P', 'N', 'R', 'R', m, n, &
717 x, ldx, work, z, ldz, w, ldw, &
718 numrnk, iwork, -1, rdummy, &
719 -1, rdummy2, -1, info1 )
720 iminwr = iwork(1)
721 mwrsvq = int(rdummy(2))
722 mlwork = max(mlwork,n+mwrsvq+int(rdummy2(1)))
723 IF ( lquery ) THEN
724 lwrsvq = int(rdummy(1))
725 olwork = max(olwork,n+lwrsvq+int(rdummy2(1)))
726 END IF
727 CASE (4)
728 jsvopt = 'J'
729 !MWRSVJ = MAX( 7, 2*M+N, 6*N+2*N*N )! for JSVOPT='V'
730 mwrsvj = max( 7, 2*m+n, 4*n+n*n, 2*n+n*n+6 )
731 mlwork = max(mlwork,n+mwrsvj)
732 iminwr = max( 3, m+3*n )
733 IF ( lquery ) THEN
734 olwork = max(olwork,n+mwrsvj)
735 END IF
736 END SELECT
737 IF ( wntvec .OR. wntex .OR. lsame(jobz,'F') ) THEN
738 jobzl = 'V'
739 ELSE
740 jobzl = 'N'
741 END IF
742 ! Workspace calculation to the SGEEV call
743 IF ( lsame(jobzl,'V') ) THEN
744 mwrkev = max( 1, 4*n )
745 ELSE
746 mwrkev = max( 1, 3*n )
747 END IF
748 mlwork = max(mlwork,n+mwrkev)
749 IF ( lquery ) THEN
750 CALL sgeev( 'N', jobzl, n, s, lds, reig, &
751 imeig, w, ldw, w, ldw, rdummy, -1, info1 )
752 lwrkev = max( mwrkev, int(rdummy(1)) )
753 olwork = max( olwork, n+lwrkev )
754 END IF
755!
756 IF ( liwork < iminwr .AND. (.NOT.lquery) ) info = -29
757 IF ( lwork < mlwork .AND. (.NOT.lquery) ) info = -27
758 END IF
759!
760 IF( info /= 0 ) THEN
761 CALL xerbla( 'SGEDMD', -info )
762 RETURN
763 ELSE IF ( lquery ) THEN
764! Return minimal and optimal workspace sizes
765 iwork(1) = iminwr
766 work(1) = mlwork
767 work(2) = olwork
768 RETURN
769 END IF
770!............................................................
771!
772 ofl = slamch('O')
773 small = slamch('S')
774 badxy = .false.
775!
776! <1> Optional scaling of the snapshots (columns of X, Y)
777! ==========================================================
778 IF ( sccolx ) THEN
779 ! The columns of X will be normalized.
780 ! To prevent overflows, the column norms of X are
781 ! carefully computed using SLASSQ.
782 k = 0
783 DO i = 1, n
784 !WORK(i) = DNRM2( M, X(1,i), 1 )
785 scale = zero
786 CALL slassq( m, x(1,i), 1, scale, ssum )
787 IF ( sisnan(scale) .OR. sisnan(ssum) ) THEN
788 k = 0
789 info = -8
790 CALL xerbla('SGEDMD',-info)
791 END IF
792 IF ( (scale /= zero) .AND. (ssum /= zero) ) THEN
793 rootsc = sqrt(ssum)
794 IF ( scale .GE. (ofl / rootsc) ) THEN
795! Norm of X(:,i) overflows. First, X(:,i)
796! is scaled by
797! ( ONE / ROOTSC ) / SCALE = 1/||X(:,i)||_2.
798! Next, the norm of X(:,i) is stored without
799! overflow as WORK(i) = - SCALE * (ROOTSC/M),
800! the minus sign indicating the 1/M factor.
801! Scaling is performed without overflow, and
802! underflow may occur in the smallest entries
803! of X(:,i). The relative backward and forward
804! errors are small in the ell_2 norm.
805 CALL slascl( 'G', 0, 0, scale, one/rootsc, &
806 m, 1, x(1,i), m, info2 )
807 work(i) = - scale * ( rootsc / float(m) )
808 ELSE
809! X(:,i) will be scaled to unit 2-norm
810 work(i) = scale * rootsc
811 CALL slascl( 'G',0, 0, work(i), one, m, 1, &
812 x(1,i), m, info2 ) ! LAPACK CALL
813! X(1:M,i) = (ONE/WORK(i)) * X(1:M,i) ! INTRINSIC
814 END IF
815 ELSE
816 work(i) = zero
817 k = k + 1
818 END IF
819 END DO
820 IF ( k == n ) THEN
821 ! All columns of X are zero. Return error code -8.
822 ! (the 8th input variable had an illegal value)
823 k = 0
824 info = -8
825 CALL xerbla('SGEDMD',-info)
826 RETURN
827 END IF
828 DO i = 1, n
829! Now, apply the same scaling to the columns of Y.
830 IF ( work(i) > zero ) THEN
831 CALL sscal( m, one/work(i), y(1,i), 1 ) ! BLAS CALL
832! Y(1:M,i) = (ONE/WORK(i)) * Y(1:M,i) ! INTRINSIC
833 ELSE IF ( work(i) < zero ) THEN
834 CALL slascl( 'G', 0, 0, -work(i), &
835 one/float(m), m, 1, y(1,i), m, info2 ) ! LAPACK CALL
836 ELSE IF ( y(isamax(m, y(1,i),1),i ) &
837 /= zero ) THEN
838! X(:,i) is zero vector. For consistency,
839! Y(:,i) should also be zero. If Y(:,i) is not
840! zero, then the data might be inconsistent or
841! corrupted. If JOBS == 'C', Y(:,i) is set to
842! zero and a warning flag is raised.
843! The computation continues but the
844! situation will be reported in the output.
845 badxy = .true.
846 IF ( lsame(jobs,'C')) &
847 CALL sscal( m, zero, y(1,i), 1 ) ! BLAS CALL
848 END IF
849 END DO
850 END IF
851 !
852 IF ( sccoly ) THEN
853 ! The columns of Y will be normalized.
854 ! To prevent overflows, the column norms of Y are
855 ! carefully computed using SLASSQ.
856 DO i = 1, n
857 !WORK(i) = DNRM2( M, Y(1,i), 1 )
858 scale = zero
859 CALL slassq( m, y(1,i), 1, scale, ssum )
860 IF ( sisnan(scale) .OR. sisnan(ssum) ) THEN
861 k = 0
862 info = -10
863 CALL xerbla('SGEDMD',-info)
864 END IF
865 IF ( scale /= zero .AND. (ssum /= zero) ) THEN
866 rootsc = sqrt(ssum)
867 IF ( scale .GE. (ofl / rootsc) ) THEN
868! Norm of Y(:,i) overflows. First, Y(:,i)
869! is scaled by
870! ( ONE / ROOTSC ) / SCALE = 1/||Y(:,i)||_2.
871! Next, the norm of Y(:,i) is stored without
872! overflow as WORK(i) = - SCALE * (ROOTSC/M),
873! the minus sign indicating the 1/M factor.
874! Scaling is performed without overflow, and
875! underflow may occur in the smallest entries
876! of Y(:,i). The relative backward and forward
877! errors are small in the ell_2 norm.
878 CALL slascl( 'G', 0, 0, scale, one/rootsc, &
879 m, 1, y(1,i), m, info2 )
880 work(i) = - scale * ( rootsc / float(m) )
881 ELSE
882! X(:,i) will be scaled to unit 2-norm
883 work(i) = scale * rootsc
884 CALL slascl( 'G',0, 0, work(i), one, m, 1, &
885 y(1,i), m, info2 ) ! LAPACK CALL
886! Y(1:M,i) = (ONE/WORK(i)) * Y(1:M,i) ! INTRINSIC
887 END IF
888 ELSE
889 work(i) = zero
890 END IF
891 END DO
892 DO i = 1, n
893! Now, apply the same scaling to the columns of X.
894 IF ( work(i) > zero ) THEN
895 CALL sscal( m, one/work(i), x(1,i), 1 ) ! BLAS CALL
896! X(1:M,i) = (ONE/WORK(i)) * X(1:M,i) ! INTRINSIC
897 ELSE IF ( work(i) < zero ) THEN
898 CALL slascl( 'G', 0, 0, -work(i), &
899 one/float(m), m, 1, x(1,i), m, info2 ) ! LAPACK CALL
900 ELSE IF ( x(isamax(m, x(1,i),1),i ) &
901 /= zero ) THEN
902! Y(:,i) is zero vector. If X(:,i) is not
903! zero, then a warning flag is raised.
904! The computation continues but the
905! situation will be reported in the output.
906 badxy = .true.
907 END IF
908 END DO
909 END IF
910!
911! <2> SVD of the data snapshot matrix X.
912! =====================================
913! The left singular vectors are stored in the array X.
914! The right singular vectors are in the array W.
915! The array W will later on contain the eigenvectors
916! of a Rayleigh quotient.
917 numrnk = n
918 SELECT CASE ( whtsvd )
919 CASE (1)
920 CALL sgesvd( 'O', 'S', m, n, x, ldx, work, b, &
921 ldb, w, ldw, work(n+1), lwork-n, info1 ) ! LAPACK CALL
922 t_or_n = 'T'
923 CASE (2)
924 CALL sgesdd( 'O', m, n, x, ldx, work, b, ldb, w, &
925 ldw, work(n+1), lwork-n, iwork, info1 ) ! LAPACK CALL
926 t_or_n = 'T'
927 CASE (3)
928 CALL sgesvdq( 'H', 'P', 'N', 'R', 'R', m, n, &
929 x, ldx, work, z, ldz, w, ldw, &
930 numrnk, iwork, liwork, work(n+max(2,m)+1),&
931 lwork-n-max(2,m), work(n+1), max(2,m), info1) ! LAPACK CALL
932 CALL slacpy( 'A', m, numrnk, z, ldz, x, ldx ) ! LAPACK CALL
933 t_or_n = 'T'
934 CASE (4)
935 CALL sgejsv( 'F', 'U', jsvopt, 'N', 'N', 'P', m, &
936 n, x, ldx, work, z, ldz, w, ldw, &
937 work(n+1), lwork-n, iwork, info1 ) ! LAPACK CALL
938 CALL slacpy( 'A', m, n, z, ldz, x, ldx ) ! LAPACK CALL
939 t_or_n = 'N'
940 xscl1 = work(n+1)
941 xscl2 = work(n+2)
942 IF ( xscl1 /= xscl2 ) THEN
943 ! This is an exceptional situation. If the
944 ! data matrices are not scaled and the
945 ! largest singular value of X overflows.
946 ! In that case SGEJSV can return the SVD
947 ! in scaled form. The scaling factor can be used
948 ! to rescale the data (X and Y).
949 CALL slascl( 'G', 0, 0, xscl1, xscl2, m, n, y, ldy, info2 )
950 END IF
951 END SELECT
952!
953 IF ( info1 > 0 ) THEN
954 ! The SVD selected subroutine did not converge.
955 ! Return with an error code.
956 info = 2
957 RETURN
958 END IF
959!
960 IF ( work(1) == zero ) THEN
961 ! The largest computed singular value of (scaled)
962 ! X is zero. Return error code -8
963 ! (the 8th input variable had an illegal value).
964 k = 0
965 info = -8
966 CALL xerbla('SGEDMD',-info)
967 RETURN
968 END IF
969!
970
971 ! snapshots matrix X. This depends on the
972 ! parameters NRNK and TOL.
973
974 SELECT CASE ( nrnk )
975 CASE ( -1 )
976 k = 1
977 DO i = 2, numrnk
978 IF ( ( work(i) <= work(1)*tol ) .OR. &
979 ( work(i) <= small ) ) EXIT
980 k = k + 1
981 END DO
982 CASE ( -2 )
983 k = 1
984 DO i = 1, numrnk-1
985 IF ( ( work(i+1) <= work(i)*tol ) .OR. &
986 ( work(i) <= small ) ) EXIT
987 k = k + 1
988 END DO
989 CASE DEFAULT
990 k = 1
991 DO i = 2, nrnk
992 IF ( work(i) <= small ) EXIT
993 k = k + 1
994 END DO
995 END SELECT
996 ! Now, U = X(1:M,1:K) is the SVD/POD basis for the
997 ! snapshot data in the input matrix X.
998
999
1000 ! Depending on the requested outputs, the computation
1001 ! is organized to compute additional auxiliary
1002 ! matrices (for the residuals and refinements).
1003 !
1004 ! In all formulas below, we need V_k*Sigma_k^(-1)
1005 ! where either V_k is in W(1:N,1:K), or V_k^T is in
1006 ! W(1:K,1:N). Here Sigma_k=diag(WORK(1:K)).
1007 IF ( lsame(t_or_n, 'N') ) THEN
1008 DO i = 1, k
1009 CALL sscal( n, one/work(i), w(1,i), 1 ) ! BLAS CALL
1010 ! W(1:N,i) = (ONE/WORK(i)) * W(1:N,i) ! INTRINSIC
1011 END DO
1012 ELSE
1013 ! This non-unit stride access is due to the fact
1014 ! that SGESVD, SGESVDQ and SGESDD return the
1015 ! transposed matrix of the right singular vectors.
1016 !DO i = 1, K
1017 ! CALL SSCAL( N, ONE/WORK(i), W(i,1), LDW ) ! BLAS CALL
1018 ! ! W(i,1:N) = (ONE/WORK(i)) * W(i,1:N) ! INTRINSIC
1019 !END DO
1020 DO i = 1, k
1021 work(n+i) = one/work(i)
1022 END DO
1023 DO j = 1, n
1024 DO i = 1, k
1025 w(i,j) = (work(n+i))*w(i,j)
1026 END DO
1027 END DO
1028 END IF
1029!
1030 IF ( wntref ) THEN
1031 !
1032 ! Need A*U(:,1:K)=Y*V_k*inv(diag(WORK(1:K)))
1033 ! for computing the refined Ritz vectors
1034 ! (optionally, outside SGEDMD).
1035 CALL sgemm( 'N', t_or_n, m, k, n, one, y, ldy, w, &
1036 ldw, zero, z, ldz ) ! BLAS CALL
1037 ! Z(1:M,1:K)=MATMUL(Y(1:M,1:N),TRANSPOSE(W(1:K,1:N))) ! INTRINSIC, for T_OR_N=='T'
1038 ! Z(1:M,1:K)=MATMUL(Y(1:M,1:N),W(1:N,1:K)) ! INTRINSIC, for T_OR_N=='N'
1039 !
1040 ! At this point Z contains
1041 ! A * U(:,1:K) = Y * V_k * Sigma_k^(-1), and
1042 ! this is needed for computing the residuals.
1043 ! This matrix is returned in the array B and
1044 ! it can be used to compute refined Ritz vectors.
1045 CALL slacpy( 'A', m, k, z, ldz, b, ldb ) ! BLAS CALL
1046 ! B(1:M,1:K) = Z(1:M,1:K) ! INTRINSIC
1047
1048 CALL sgemm( 'T', 'N', k, k, m, one, x, ldx, z, &
1049 ldz, zero, s, lds ) ! BLAS CALL
1050 ! S(1:K,1:K) = MATMUL(TANSPOSE(X(1:M,1:K)),Z(1:M,1:K)) ! INTRINSIC
1051 ! At this point S = U^T * A * U is the Rayleigh quotient.
1052 ELSE
1053 ! A * U(:,1:K) is not explicitly needed and the
1054 ! computation is organized differently. The Rayleigh
1055 ! quotient is computed more efficiently.
1056 CALL sgemm( 'T', 'N', k, n, m, one, x, ldx, y, ldy, &
1057 zero, z, ldz ) ! BLAS CALL
1058 ! Z(1:K,1:N) = MATMUL( TRANSPOSE(X(1:M,1:K)), Y(1:M,1:N) ) ! INTRINSIC
1059 ! In the two SGEMM calls here, can use K for LDZ
1060 CALL sgemm( 'N', t_or_n, k, k, n, one, z, ldz, w, &
1061 ldw, zero, s, lds ) ! BLAS CALL
1062 ! S(1:K,1:K) = MATMUL(Z(1:K,1:N),TRANSPOSE(W(1:K,1:N))) ! INTRINSIC, for T_OR_N=='T'
1063 ! S(1:K,1:K) = MATMUL(Z(1:K,1:N),(W(1:N,1:K))) ! INTRINSIC, for T_OR_N=='N'
1064 ! At this point S = U^T * A * U is the Rayleigh quotient.
1065 ! If the residuals are requested, save scaled V_k into Z.
1066 ! Recall that V_k or V_k^T is stored in W.
1067 IF ( wntres .OR. wntex ) THEN
1068 IF ( lsame(t_or_n, 'N') ) THEN
1069 CALL slacpy( 'A', n, k, w, ldw, z, ldz )
1070 ELSE
1071 CALL slacpy( 'A', k, n, w, ldw, z, ldz )
1072 END IF
1073 END IF
1074 END IF
1075!
1076
1077 ! right eigenvectors of the Rayleigh quotient.
1078 !
1079 CALL sgeev( 'N', jobzl, k, s, lds, reig, imeig, w, &
1080 ldw, w, ldw, work(n+1), lwork-n, info1 ) ! LAPACK CALL
1081 !
1082 ! W(1:K,1:K) contains the eigenvectors of the Rayleigh
1083 ! quotient. Even in the case of complex spectrum, all
1084 ! computation is done in real arithmetic. REIG and
1085 ! IMEIG are the real and the imaginary parts of the
1086 ! eigenvalues, so that the spectrum is given as
1087 ! REIG(:) + sqrt(-1)*IMEIG(:). Complex conjugate pairs
1088 ! are listed at consecutive positions. For such a
1089 ! complex conjugate pair of the eigenvalues, the
1090 ! corresponding eigenvectors are also a complex
1091 ! conjugate pair with the real and imaginary parts
1092 ! stored column-wise in W at the corresponding
1093 ! consecutive column indices. See the description of Z.
1094 ! Also, see the description of SGEEV.
1095 IF ( info1 > 0 ) THEN
1096 ! SGEEV failed to compute the eigenvalues and
1097 ! eigenvectors of the Rayleigh quotient.
1098 info = 3
1099 RETURN
1100 END IF
1101!
1102 ! <6> Compute the eigenvectors (if requested) and,
1103 ! the residuals (if requested).
1104 !
1105 IF ( wntvec .OR. wntex ) THEN
1106 IF ( wntres ) THEN
1107 IF ( wntref ) THEN
1108 ! Here, if the refinement is requested, we have
1109 ! A*U(:,1:K) already computed and stored in Z.
1110 ! For the residuals, need Y = A * U(:,1;K) * W.
1111 CALL sgemm( 'N', 'N', m, k, k, one, z, ldz, w, &
1112 ldw, zero, y, ldy ) ! BLAS CALL
1113 ! Y(1:M,1:K) = Z(1:M,1:K) * W(1:K,1:K) ! INTRINSIC
1114 ! This frees Z; Y contains A * U(:,1:K) * W.
1115 ELSE
1116 ! Compute S = V_k * Sigma_k^(-1) * W, where
1117 ! V_k * Sigma_k^(-1) is stored in Z
1118 CALL sgemm( t_or_n, 'N', n, k, k, one, z, ldz, &
1119 w, ldw, zero, s, lds )
1120 ! Then, compute Z = Y * S =
1121 ! = Y * V_k * Sigma_k^(-1) * W(1:K,1:K) =
1122 ! = A * U(:,1:K) * W(1:K,1:K)
1123 CALL sgemm( 'N', 'N', m, k, n, one, y, ldy, s, &
1124 lds, zero, z, ldz )
1125 ! Save a copy of Z into Y and free Z for holding
1126 ! the Ritz vectors.
1127 CALL slacpy( 'A', m, k, z, ldz, y, ldy )
1128 IF ( wntex ) CALL slacpy( 'A', m, k, z, ldz, b, ldb )
1129 END IF
1130 ELSE IF ( wntex ) THEN
1131 ! Compute S = V_k * Sigma_k^(-1) * W, where
1132 ! V_k * Sigma_k^(-1) is stored in Z
1133 CALL sgemm( t_or_n, 'N', n, k, k, one, z, ldz, &
1134 w, ldw, zero, s, lds )
1135 ! Then, compute Z = Y * S =
1136 ! = Y * V_k * Sigma_k^(-1) * W(1:K,1:K) =
1137 ! = A * U(:,1:K) * W(1:K,1:K)
1138 CALL sgemm( 'N', 'N', m, k, n, one, y, ldy, s, &
1139 lds, zero, b, ldb )
1140 ! The above call replaces the following two calls
1141 ! that were used in the developing-testing phase.
1142 ! CALL SGEMM( 'N', 'N', M, K, N, ONE, Y, LDY, S, &
1143 ! LDS, ZERO, Z, LDZ)
1144 ! Save a copy of Z into B and free Z for holding
1145 ! the Ritz vectors.
1146 ! CALL SLACPY( 'A', M, K, Z, LDZ, B, LDB )
1147 END IF
1148!
1149 ! Compute the real form of the Ritz vectors
1150 IF ( wntvec ) CALL sgemm( 'N', 'N', m, k, k, one, x, ldx, w, ldw, &
1151 zero, z, ldz ) ! BLAS CALL
1152 ! Z(1:M,1:K) = MATMUL(X(1:M,1:K), W(1:K,1:K)) ! INTRINSIC
1153!
1154 IF ( wntres ) THEN
1155 i = 1
1156 DO WHILE ( i <= k )
1157 IF ( imeig(i) == zero ) THEN
1158 ! have a real eigenvalue with real eigenvector
1159 CALL saxpy( m, -reig(i), z(1,i), 1, y(1,i), 1 ) ! BLAS CALL
1160 ! Y(1:M,i) = Y(1:M,i) - REIG(i) * Z(1:M,i) ! INTRINSIC
1161 res(i) = snrm2( m, y(1,i), 1 ) ! BLAS CALL
1162 i = i + 1
1163 ELSE
1164 ! Have a complex conjugate pair
1165 ! REIG(i) +- sqrt(-1)*IMEIG(i).
1166 ! Since all computation is done in real
1167 ! arithmetic, the formula for the residual
1168 ! is recast for real representation of the
1169 ! complex conjugate eigenpair. See the
1170 ! description of RES.
1171 ab(1,1) = reig(i)
1172 ab(2,1) = -imeig(i)
1173 ab(1,2) = imeig(i)
1174 ab(2,2) = reig(i)
1175 CALL sgemm( 'N', 'N', m, 2, 2, -one, z(1,i), &
1176 ldz, ab, 2, one, y(1,i), ldy ) ! BLAS CALL
1177 ! Y(1:M,i:i+1) = Y(1:M,i:i+1) - Z(1:M,i:i+1) * AB ! INTRINSIC
1178 res(i) = slange( 'F', m, 2, y(1,i), ldy, &
1179 work(n+1) ) ! LAPACK CALL
1180 res(i+1) = res(i)
1181 i = i + 2
1182 END IF
1183 END DO
1184 END IF
1185 END IF
1186!
1187 IF ( whtsvd == 4 ) THEN
1188 work(n+1) = xscl1
1189 work(n+2) = xscl2
1190 END IF
1191!
1192! Successful exit.
1193 IF ( .NOT. badxy ) THEN
1194 info = 0
1195 ELSE
1196 ! A warning on possible data inconsistency.
1197 ! This should be a rare event.
1198 info = 4
1199 END IF
1200!............................................................
1201 RETURN
1202! ......
1203 END SUBROUTINE sgedmd
1204
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine saxpy(n, sa, sx, incx, sy, incy)
SAXPY
Definition saxpy.f:89
subroutine sgedmd(jobs, jobz, jobr, jobf, whtsvd, m, n, x, ldx, y, ldy, nrnk, tol, k, reig, imeig, z, ldz, res, b, ldb, w, ldw, s, lds, work, lwork, iwork, liwork, info)
SGEDMD computes the Dynamic Mode Decomposition (DMD) for a pair of data snapshot matrices.
Definition sgedmd.f90:535
subroutine sgeev(jobvl, jobvr, n, a, lda, wr, wi, vl, ldvl, vr, ldvr, work, lwork, info)
SGEEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for GE matrices
Definition sgeev.f:192
subroutine sgejsv(joba, jobu, jobv, jobr, jobt, jobp, m, n, a, lda, sva, u, ldu, v, ldv, work, lwork, iwork, info)
SGEJSV
Definition sgejsv.f:476
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:188
subroutine sgesdd(jobz, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, iwork, info)
SGESDD
Definition sgesdd.f:213
subroutine sgesvd(jobu, jobvt, m, n, a, lda, s, u, ldu, vt, ldvt, work, lwork, info)
SGESVD computes the singular value decomposition (SVD) for GE matrices
Definition sgesvd.f:211
subroutine sgesvdq(joba, jobp, jobr, jobu, jobv, m, n, a, lda, s, u, ldu, v, ldv, numrank, iwork, liwork, work, lwork, rwork, lrwork, info)
SGESVDQ computes the singular value decomposition (SVD) with a QR-Preconditioned QR SVD Method for GE...
Definition sgesvdq.f:415
integer function isamax(n, sx, incx)
ISAMAX
Definition isamax.f:71
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
real function slange(norm, m, n, a, lda, work)
SLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition slange.f:114
subroutine slascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
SLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition slascl.f:143
subroutine slassq(n, x, incx, scale, sumsq)
SLASSQ updates a sum of squares represented in scaled form.
Definition slassq.f90:124
real(wp) function snrm2(n, x, incx)
SNRM2
Definition snrm2.f90:89
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79