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

◆ dgedmd()

subroutine dgedmd ( 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 
)

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

Purpose:
    DGEDMD 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, DGEDMD 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, DGEDMD 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:
    Approved for Public Release, Distribution Unlimited.
    Cleared by DARPA on September 29, 2022
Parameters
[in]JOBS
    JOBS (input) is 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 :: DGESVD (the QR SVD algorithm)
    2 :: DGESDD (the Divide and Conquer algorithm; if enough
         workspace available, this is the fastest option)
    3 :: DGESVDQ (the preconditioned QR SVD  ; this and 4
         are the most accurate options)
    4 :: DGEJSV (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
    right 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 DGEEV.
    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 DGEDMD is only workspace query, then
    WORK(1) contains the minimal workspace length and
    WORK(2) is the optimal workspace length. Hence, the
    leng 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 DGESVD.
    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 DGESDD.
    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 DGESVDQ.
    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 DGEJSV.
    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 531 of file dgedmd.f90.

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