LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ sgedmd()

subroutine sgedmd ( character, intent(in)  jobs,
character, intent(in)  jobz,
character, intent(in)  jobr,
character, intent(in)  jobf,
integer, intent(in)  whtsvd,
integer, intent(in)  m,
integer, intent(in)  n,
real(kind=wp), dimension(ldx,*), intent(inout)  x,
integer, intent(in)  ldx,
real(kind=wp), dimension(ldy,*), intent(inout)  y,
integer, intent(in)  ldy,
integer, intent(in)  nrnk,
real(kind=wp), intent(in)  tol,
integer, intent(out)  k,
real(kind=wp), dimension(*), intent(out)  reig,
real(kind=wp), dimension(*), intent(out)  imeig,
real(kind=wp), dimension(ldz,*), intent(out)  z,
integer, intent(in)  ldz,
real(kind=wp), dimension(*), intent(out)  res,
real(kind=wp), dimension(ldb,*), intent(out)  b,
integer, intent(in)  ldb,
real(kind=wp), dimension(ldw,*), intent(out)  w,
integer, intent(in)  ldw,
real(kind=wp), dimension(lds,*), intent(out)  s,
integer, intent(in)  lds,
real(kind=wp), dimension(*), intent(out)  work,
integer, intent(in)  lwork,
integer, dimension(*), intent(out)  iwork,
integer, intent(in)  liwork,
integer, intent(out)  info 
)

SGEDMD computes the Dynamic Mode Decomposition (DMD) for a pair of data snapshot matrices.

Purpose:
    SGEDMD computes the Dynamic Mode Decomposition (DMD) for
    a pair of data snapshot matrices. For the input matrices
    X and Y such that Y = A*X with an unaccessible matrix
    A, SGEDMD computes a certain number of Ritz pairs of A using
    the standard Rayleigh-Ritz extraction from a subspace of
    range(X) that is determined using the leading left singular
    vectors of X. Optionally, SGEDMD returns the residuals
    of the computed Ritz pairs, the information needed for
    a refinement of the Ritz vectors, or the eigenvectors of
    the Exact DMD.
    For further details see the references listed
    below. For more details of the implementation see [3].
References:
    [1] P. Schmid: Dynamic mode decomposition of numerical
        and experimental data,
        Journal of Fluid Mechanics 656, 5-28, 2010.
    [2] Z. Drmac, I. Mezic, R. Mohr: Data driven modal
        decompositions: analysis and enhancements,
        SIAM J. on Sci. Comp. 40 (4), A2253-A2285, 2018.
    [3] Z. Drmac: A LAPACK implementation of the Dynamic
        Mode Decomposition I. Technical report. AIMDyn Inc.
        and LAPACK Working Note 298.
    [4] J. Tu, C. W. Rowley, D. M. Luchtenburg, S. L.
        Brunton, N. Kutz: On Dynamic Mode Decomposition:
        Theory and Applications, Journal of Computational
        Dynamics 1(2), 391 -421, 2014.
Developed and supported by:
    Developed and coded by Zlatko Drmac, Faculty of Science,
    University of Zagreb;  drmac@math.hr
    In cooperation with
    AIMdyn Inc., Santa Barbara, CA.
    and supported by
    - DARPA SBIR project "Koopman Operator-Based Forecasting
    for Nonstationary Processes from Near-Term, Limited
    Observational Data" Contract No: W31P4Q-21-C-0007
    - DARPA PAI project "Physics-Informed Machine Learning
    Methodologies" Contract No: HR0011-18-9-0033
    - DARPA MoDyL project "A Data-Driven, Operator-Theoretic
    Framework for Space-Time Analysis of Process Dynamics"
    Contract No: HR0011-16-C-0116
    Any opinions, findings and conclusions or recommendations
    expressed in this material are those of the author and
    do not necessarily reflect the views of the DARPA SBIR
    Program Office
Distribution Statement A:
    Distribution Statement A:
    Approved for Public Release, Distribution Unlimited.
    Cleared by DARPA on September 29, 2022
Parameters
[in]JOBS
    JOBS (input) CHARACTER*1
    Determines whether the initial data snapshots are scaled
    by a diagonal matrix.
    'S' :: The data snapshots matrices X and Y are multiplied
           with a diagonal matrix D so that X*D has unit
           nonzero columns (in the Euclidean 2-norm)
    'C' :: The snapshots are scaled as with the 'S' option.
           If it is found that an i-th column of X is zero
           vector and the corresponding i-th column of Y is
           non-zero, then the i-th column of Y is set to
           zero and a warning flag is raised.
    'Y' :: The data snapshots matrices X and Y are multiplied
           by a diagonal matrix D so that Y*D has unit
           nonzero columns (in the Euclidean 2-norm)
    'N' :: No data scaling.
[in]JOBZ
    JOBZ (input) CHARACTER*1
    Determines whether the eigenvectors (Koopman modes) will
    be computed.
    'V' :: The eigenvectors (Koopman modes) will be computed
           and returned in the matrix Z.
           See the description of Z.
    'F' :: The eigenvectors (Koopman modes) will be returned
           in factored form as the product X(:,1:K)*W, where X
           contains a POD basis (leading left singular vectors
           of the data matrix X) and W contains the eigenvectors
           of the corresponding Rayleigh quotient.
           See the descriptions of K, X, W, Z.
    'N' :: The eigenvectors are not computed.
[in]JOBR
    JOBR (input) CHARACTER*1
    Determines whether to compute the residuals.
    'R' :: The residuals for the computed eigenpairs will be
           computed and stored in the array RES.
           See the description of RES.
           For this option to be legal, JOBZ must be 'V'.
    'N' :: The residuals are not computed.
[in]JOBF
    JOBF (input) CHARACTER*1
    Specifies whether to store information needed for post-
    processing (e.g. computing refined Ritz vectors)
    'R' :: The matrix needed for the refinement of the Ritz
           vectors is computed and stored in the array B.
           See the description of B.
    'E' :: The unscaled eigenvectors of the Exact DMD are
           computed and returned in the array B. See the
           description of B.
    'N' :: No eigenvector refinement data is computed.
[in]WHTSVD
    WHTSVD (input) INTEGER, WHSTVD in { 1, 2, 3, 4 }
    Allows for a selection of the SVD algorithm from the
    LAPACK library.
    1 :: SGESVD (the QR SVD algorithm)
    2 :: SGESDD (the Divide and Conquer algorithm; if enough
         workspace available, this is the fastest option)
    3 :: SGESVDQ (the preconditioned QR SVD  ; this and 4
         are the most accurate options)
    4 :: SGEJSV (the preconditioned Jacobi SVD; this and 3
         are the most accurate options)
    For the four methods above, a significant difference in
    the accuracy of small singular values is possible if
    the snapshots vary in norm so that X is severely
    ill-conditioned. If small (smaller than EPS*||X||)
    singular values are of interest and JOBS=='N',  then
    the options (3, 4) give the most accurate results, where
    the option 4 is slightly better and with stronger
    theoretical background.
    If JOBS=='S', i.e. the columns of X will be normalized,
    then all methods give nearly equally accurate results.
[in]M
    M (input) INTEGER, M>= 0
    The state space dimension (the row dimension of X, Y).
[in]N
    N (input) INTEGER, 0 <= N <= M
    The number of data snapshot pairs
    (the number of columns of X and Y).
[in,out]X
    X (input/output) REAL(KIND=WP) M-by-N array
    > On entry, X contains the data snapshot matrix X. It is
    assumed that the column norms of X are in the range of
    the normalized floating point numbers.
    < On exit, the leading K columns of X contain a POD basis,
    i.e. the leading K left singular vectors of the input
    data matrix X, U(:,1:K). All N columns of X contain all
    left singular vectors of the input matrix X.
    See the descriptions of K, Z and W.
[in]LDX
    LDX (input) INTEGER, LDX >= M
    The leading dimension of the array X.
[in,out]Y
    Y (input/workspace/output) REAL(KIND=WP) M-by-N array
    > On entry, Y contains the data snapshot matrix Y
    < On exit,
    If JOBR == 'R', the leading K columns of Y  contain
    the residual vectors for the computed Ritz pairs.
    See the description of RES.
    If JOBR == 'N', Y contains the original input data,
                    scaled according to the value of JOBS.
[in]LDY
    LDY (input) INTEGER , LDY >= M
    The leading dimension of the array Y.
[in]NRNK
    NRNK (input) INTEGER
    Determines the mode how to compute the numerical rank,
    i.e. how to truncate small singular values of the input
    matrix X. On input, if
    NRNK = -1 :: i-th singular value sigma(i) is truncated
                 if sigma(i) <= TOL*sigma(1)
                 This option is recommended.
    NRNK = -2 :: i-th singular value sigma(i) is truncated
                 if sigma(i) <= TOL*sigma(i-1)
                 This option is included for R&D purposes.
                 It requires highly accurate SVD, which
                 may not be feasible.
    The numerical rank can be enforced by using positive
    value of NRNK as follows:
    0 < NRNK <= N :: at most NRNK largest singular values
    will be used. If the number of the computed nonzero
    singular values is less than NRNK, then only those
    nonzero values will be used and the actually used
    dimension is less than NRNK. The actual number of
    the nonzero singular values is returned in the variable
    K. See the descriptions of TOL and  K.
[in]TOL
    TOL (input) REAL(KIND=WP), 0 <= TOL < 1
    The tolerance for truncating small singular values.
    See the description of NRNK.
[out]K
    K (output) INTEGER,  0 <= K <= N
    The dimension of the POD basis for the data snapshot
    matrix X and the number of the computed Ritz pairs.
    The value of K is determined according to the rule set
    by the parameters NRNK and TOL.
    See the descriptions of NRNK and TOL.
[out]REIG
    REIG (output) REAL(KIND=WP) N-by-1 array
    The leading K (K<=N) entries of REIG contain
    the real parts of the computed eigenvalues
    REIG(1:K) + sqrt(-1)*IMEIG(1:K).
    See the descriptions of K, IMEIG, and Z.
[out]IMEIG
    IMEIG (output) REAL(KIND=WP) N-by-1 array
    The leading K (K<=N) entries of IMEIG contain
    the imaginary parts of the computed eigenvalues
    REIG(1:K) + sqrt(-1)*IMEIG(1:K).
    The eigenvalues are determined as follows:
    If IMEIG(i) == 0, then the corresponding eigenvalue is
    real, LAMBDA(i) = REIG(i).
    If IMEIG(i)>0, then the corresponding complex
    conjugate pair of eigenvalues reads
    LAMBDA(i)   = REIG(i) + sqrt(-1)*IMAG(i)
    LAMBDA(i+1) = REIG(i) - sqrt(-1)*IMAG(i)
    That is, complex conjugate pairs have consecutive
    indices (i,i+1), with the positive imaginary part
    listed first.
    See the descriptions of K, REIG, and Z.
[out]Z
    Z (workspace/output) REAL(KIND=WP)  M-by-N array
    If JOBZ =='V' then
       Z contains real Ritz vectors as follows:
       If IMEIG(i)=0, then Z(:,i) is an eigenvector of
       the i-th Ritz value; ||Z(:,i)||_2=1.
       If IMEIG(i) > 0 (and IMEIG(i+1) < 0) then
       [Z(:,i) Z(:,i+1)] span an invariant subspace and
       the Ritz values extracted from this subspace are
       REIG(i) + sqrt(-1)*IMEIG(i) and
       REIG(i) - sqrt(-1)*IMEIG(i).
       The corresponding eigenvectors are
       Z(:,i) + sqrt(-1)*Z(:,i+1) and
       Z(:,i) - sqrt(-1)*Z(:,i+1), respectively.
       || Z(:,i:i+1)||_F = 1.
    If JOBZ == 'F', then the above descriptions hold for
    the columns of X(:,1:K)*W(1:K,1:K), where the columns
    of W(1:k,1:K) are the computed eigenvectors of the
    K-by-K Rayleigh quotient. The columns of W(1:K,1:K)
    are similarly structured: If IMEIG(i) == 0 then
    X(:,1:K)*W(:,i) is an eigenvector, and if IMEIG(i)>0
    then X(:,1:K)*W(:,i)+sqrt(-1)*X(:,1:K)*W(:,i+1) and
         X(:,1:K)*W(:,i)-sqrt(-1)*X(:,1:K)*W(:,i+1)
    are the eigenvectors of LAMBDA(i), LAMBDA(i+1).
    See the descriptions of REIG, IMEIG, X and W.
[in]LDZ
    LDZ (input) INTEGER , LDZ >= M
    The leading dimension of the array Z.
[out]RES
    RES (output) REAL(KIND=WP) N-by-1 array
    RES(1:K) contains the residuals for the K computed
    Ritz pairs.
    If LAMBDA(i) is real, then
       RES(i) = || A * Z(:,i) - LAMBDA(i)*Z(:,i))||_2.
    If [LAMBDA(i), LAMBDA(i+1)] is a complex conjugate pair
    then
    RES(i)=RES(i+1) = || A * Z(:,i:i+1) - Z(:,i:i+1) *B||_F
    where B = [ real(LAMBDA(i)) imag(LAMBDA(i)) ]
              [-imag(LAMBDA(i)) real(LAMBDA(i)) ].
    It holds that
    RES(i)   = || A*ZC(:,i)   - LAMBDA(i)  *ZC(:,i)   ||_2
    RES(i+1) = || A*ZC(:,i+1) - LAMBDA(i+1)*ZC(:,i+1) ||_2
    where ZC(:,i)   =  Z(:,i) + sqrt(-1)*Z(:,i+1)
          ZC(:,i+1) =  Z(:,i) - sqrt(-1)*Z(:,i+1)
    See the description of REIG, IMEIG and Z.
[out]B
    B (output) REAL(KIND=WP)  M-by-N array.
    IF JOBF =='R', B(1:M,1:K) contains A*U(:,1:K), and can
    be used for computing the refined vectors; see further
    details in the provided references.
    If JOBF == 'E', B(1:M,1;K) contains
    A*U(:,1:K)*W(1:K,1:K), which are the vectors from the
    Exact DMD, up to scaling by the inverse eigenvalues.
    If JOBF =='N', then B is not referenced.
    See the descriptions of X, W, K.
[in]LDB
    LDB (input) INTEGER, LDB >= M
    The leading dimension of the array B.
[out]W
    W (workspace/output) REAL(KIND=WP) N-by-N array
    On exit, W(1:K,1:K) contains the K computed
    eigenvectors of the matrix Rayleigh quotient (real and
    imaginary parts for each complex conjugate pair of the
    eigenvalues). The Ritz vectors (returned in Z) are the
    product of X (containing a POD basis for the input
    matrix X) and W. See the descriptions of K, S, X and Z.
    W is also used as a workspace to temporarily store the
    left singular vectors of X.
[in]LDW
    LDW (input) INTEGER, LDW >= N
    The leading dimension of the array W.
[out]S
    S (workspace/output) REAL(KIND=WP) N-by-N array
    The array S(1:K,1:K) is used for the matrix Rayleigh
    quotient. This content is overwritten during
    the eigenvalue decomposition by SGEEV.
    See the description of K.
[in]LDS
    LDS (input) INTEGER, LDS >= N
    The leading dimension of the array S.
[out]WORK
    WORK (workspace/output) REAL(KIND=WP) LWORK-by-1 array
    On exit, WORK(1:N) contains the singular values of
    X (for JOBS=='N') or column scaled X (JOBS=='S', 'C').
    If WHTSVD==4, then WORK(N+1) and WORK(N+2) contain
    scaling factor WORK(N+2)/WORK(N+1) used to scale X
    and Y to avoid overflow in the SVD of X.
    This may be of interest if the scaling option is off
    and as many as possible smallest eigenvalues are
    desired to the highest feasible accuracy.
    If the call to SGEDMD is only workspace query, then
    WORK(1) contains the minimal workspace length and
    WORK(2) is the optimal workspace length. Hence, the
    length of work is at least 2.
    See the description of LWORK.
[in]LWORK
    LWORK (input) INTEGER
    The minimal length of the workspace vector WORK.
    LWORK is calculated as follows:
    If WHTSVD == 1 ::
       If JOBZ == 'V', then
       LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,4*N)).
       If JOBZ == 'N'  then
       LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,3*N)).
       Here LWORK_SVD = MAX(1,3*N+M,5*N) is the minimal
       workspace length of SGESVD.
    If WHTSVD == 2 ::
       If JOBZ == 'V', then
       LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,4*N))
       If JOBZ == 'N', then
       LWORK >= MAX(2, N + LWORK_SVD, N+MAX(1,3*N))
       Here LWORK_SVD = MAX(M, 5*N*N+4*N)+3*N*N is the
       minimal workspace length of SGESDD.
    If WHTSVD == 3 ::
       If JOBZ == 'V', then
       LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,4*N))
       If JOBZ == 'N', then
       LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,3*N))
       Here LWORK_SVD = N+M+MAX(3*N+1,
                       MAX(1,3*N+M,5*N),MAX(1,N))
       is the minimal workspace length of SGESVDQ.
    If WHTSVD == 4 ::
       If JOBZ == 'V', then
       LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,4*N))
       If JOBZ == 'N', then
       LWORK >= MAX(2, N+LWORK_SVD,N+MAX(1,3*N))
       Here LWORK_SVD = MAX(7,2*M+N,6*N+2*N*N) is the
       minimal workspace length of SGEJSV.
    The above expressions are not simplified in order to
    make the usage of WORK more transparent, and for
    easier checking. In any case, LWORK >= 2.
    If on entry LWORK = -1, then a workspace query is
    assumed and the procedure only computes the minimal
    and the optimal workspace lengths for both WORK and
    IWORK. See the descriptions of WORK and IWORK.
[out]IWORK
    IWORK (workspace/output) INTEGER LIWORK-by-1 array
    Workspace that is required only if WHTSVD equals
    2 , 3 or 4. (See the description of WHTSVD).
    If on entry LWORK =-1 or LIWORK=-1, then the
    minimal length of IWORK is computed and returned in
    IWORK(1). See the description of LIWORK.
[in]LIWORK
    LIWORK (input) INTEGER
    The minimal length of the workspace vector IWORK.
    If WHTSVD == 1, then only IWORK(1) is used; LIWORK >=1
    If WHTSVD == 2, then LIWORK >= MAX(1,8*MIN(M,N))
    If WHTSVD == 3, then LIWORK >= MAX(1,M+N-1)
    If WHTSVD == 4, then LIWORK >= MAX(3,M+3*N)
    If on entry LIWORK = -1, then a workspace query is
    assumed and the procedure only computes the minimal
    and the optimal workspace lengths for both WORK and
    IWORK. See the descriptions of WORK and IWORK.
[out]INFO
    INFO (output) INTEGER
    -i < 0 :: On entry, the i-th argument had an
              illegal value
       = 0 :: Successful return.
       = 1 :: Void input. Quick exit (M=0 or N=0).
       = 2 :: The SVD computation of X did not converge.
              Suggestion: Check the input data and/or
              repeat with different WHTSVD.
       = 3 :: The computation of the eigenvalues did not
              converge.
       = 4 :: If data scaling was requested on input and
              the procedure found inconsistency in the data
              such that for some column index i,
              X(:,i) = 0 but Y(:,i) /= 0, then Y(:,i) is set
              to zero if JOBS=='C'. The computation proceeds
              with original or modified data and warning
              flag is set with INFO=4.
Author
Zlatko Drmac

Definition at line 530 of file sgedmd.f90.

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! ......
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine saxpy(n, sa, sx, incx, sy, incy)
SAXPY
Definition saxpy.f:89
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
logical function sisnan(sin)
SISNAN tests input for NaN.
Definition sisnan.f:59
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
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real(wp) function snrm2(n, x, incx)
SNRM2
Definition snrm2.f90:89
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
Here is the call graph for this function:
Here is the caller graph for this function: