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

◆ zgedmd()

subroutine zgedmd ( 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,
complex(kind=wp), dimension(ldx,*), intent(inout)  x,
integer, intent(in)  ldx,
complex(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,
complex(kind=wp), dimension(*), intent(out)  eigs,
complex(kind=wp), dimension(ldz,*), intent(out)  z,
integer, intent(in)  ldz,
real(kind=wp), dimension(*), intent(out)  res,
complex(kind=wp), dimension(ldb,*), intent(out)  b,
integer, intent(in)  ldb,
complex(kind=wp), dimension(ldw,*), intent(out)  w,
integer, intent(in)  ldw,
complex(kind=wp), dimension(lds,*), intent(out)  s,
integer, intent(in)  lds,
complex(kind=wp), dimension(*), intent(out)  zwork,
integer, intent(in)  lzwork,
real(kind=wp), dimension(*), intent(out)  rwork,
integer, intent(in)  lrwork,
integer, dimension(*), intent(out)  iwork,
integer, intent(in)  liwork,
integer, intent(out)  info 
)

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

Purpose:
    ZGEDMD 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, ZGEDMD 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, ZGEDMD 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) 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 :: ZGESVD (the QR SVD algorithm)
    2 :: ZGESDD (the Divide and Conquer algorithm; if enough
         workspace available, this is the fastest option)
    3 :: ZGESVDQ (the preconditioned QR SVD  ; this and 4
         are the most accurate options)
    4 :: ZGEJSV (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]LDX
    X (input/output) COMPLEX(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.
LDX (input) INTEGER, LDX >= M
    The leading dimension of the array X.
[in,out]Y
    Y (input/workspace/output) COMPLEX(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]EIGS
    EIGS (output) COMPLEX(KIND=WP) N-by-1 array
    The leading K (K<=N) entries of EIGS contain
    the computed eigenvalues (Ritz values).
    See the descriptions of K, and Z.
[out]Z
    Z (workspace/output) COMPLEX(KIND=WP)  M-by-N array
    If JOBZ =='V' then Z contains the  Ritz vectors.  Z(:,i)
    is an eigenvector of the i-th Ritz value; ||Z(:,i)||_2=1.
    If JOBZ == 'F', then the Z(:,i)'s are given implicitly as
    the columns of X(:,1:K)*W(1:K,1:K), i.e. X(:,1:K)*W(:,i)
    is an eigenvector corresponding to EIGS(i). The columns
    of W(1:k,1:K) are the computed eigenvectors of the
    K-by-K Rayleigh quotient.
    See the descriptions of EIGS, 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,
    RES(i) = || A * Z(:,i) - EIGS(i)*Z(:,i))||_2.
    See the description of EIGS and Z.
[out]B
    B (output) COMPLEX(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) COMPLEX(KIND=WP) N-by-N array
    On exit, W(1:K,1:K) contains the K computed
    eigenvectors of the matrix Rayleigh quotient.
    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) COMPLEX(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 ZGEEV.
    See the description of K.
[in]LDS
    LDS (input) INTEGER, LDS >= N
    The leading dimension of the array S.
[out]ZWORK
    ZWORK (workspace/output) COMPLEX(KIND=WP) LZWORK-by-1 array
    ZWORK is used as complex workspace in the complex SVD, as
    specified by WHTSVD (1,2, 3 or 4) and for ZGEEV for computing
    the eigenvalues of a Rayleigh quotient.
    If the call to ZGEDMD is only workspace query, then
    ZWORK(1) contains the minimal complex workspace length and
    ZWORK(2) is the optimal complex workspace length.
    Hence, the length of work is at least 2.
    See the description of LZWORK.
[in]LZWORK
    LZWORK (input) INTEGER
    The minimal length of the workspace vector ZWORK.
    LZWORK is calculated as MAX(LZWORK_SVD, LZWORK_ZGEEV),
    where LZWORK_ZGEEV = MAX( 1, 2*N )  and the minimal
    LZWORK_SVD is calculated as follows
    If WHTSVD == 1 :: ZGESVD ::
       LZWORK_SVD = MAX(1,2*MIN(M,N)+MAX(M,N))
    If WHTSVD == 2 :: ZGESDD ::
       LZWORK_SVD = 2*MIN(M,N)*MIN(M,N)+2*MIN(M,N)+MAX(M,N)
    If WHTSVD == 3 :: ZGESVDQ ::
       LZWORK_SVD = obtainable by a query
    If WHTSVD == 4 :: ZGEJSV ::
       LZWORK_SVD = obtainable by a query
    If on entry LZWORK = -1, then a workspace query is
    assumed and the procedure only computes the minimal
    and the optimal workspace lengths and returns them in
    LZWORK(1) and LZWORK(2), respectively.
[out]RWORK
    RWORK (workspace/output) REAL(KIND=WP) LRWORK-by-1 array
    On exit, RWORK(1:N) contains the singular values of
    X (for JOBS=='N') or column scaled X (JOBS=='S', 'C').
    If WHTSVD==4, then RWORK(N+1) and RWORK(N+2) contain
    scaling factor RWORK(N+2)/RWORK(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 ZGEDMD is only workspace query, then
    RWORK(1) contains the minimal workspace length.
    See the description of LRWORK.
[in]LRWORK
    LRWORK (input) INTEGER
    The minimal length of the workspace vector RWORK.
    LRWORK is calculated as follows:
    LRWORK = MAX(1, N+LRWORK_SVD,N+LRWORK_ZGEEV), where
    LRWORK_ZGEEV = MAX(1,2*N) and RWORK_SVD is the real workspace
    for the SVD subroutine determined by the input parameter
    WHTSVD.
    If WHTSVD == 1 :: ZGESVD ::
       LRWORK_SVD = 5*MIN(M,N)
    If WHTSVD == 2 :: ZGESDD ::
       LRWORK_SVD =  MAX(5*MIN(M,N)*MIN(M,N)+7*MIN(M,N),
       2*MAX(M,N)*MIN(M,N)+2*MIN(M,N)*MIN(M,N)+MIN(M,N) ) )
    If WHTSVD == 3 :: ZGESVDQ ::
       LRWORK_SVD = obtainable by a query
    If WHTSVD == 4 :: ZGEJSV ::
       LRWORK_SVD = obtainable by a query
    If on entry LRWORK = -1, then a workspace query is
    assumed and the procedure only computes the minimal
    real workspace length and returns it in RWORK(1).
[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  ZWORK, RWORK and
    IWORK. See the descriptions of ZWORK, RWORK 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 493 of file zgedmd.f90.

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