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

◆ cgedmd()

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

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

Purpose:
    CGEDMD 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, CGEDMD 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, CGEDMD 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 :: CGESVD (the QR SVD algorithm)
    2 :: CGESDD (the Divide and Conquer algorithm; if enough
         workspace available, this is the fastest option)
    3 :: CGESVDQ (the preconditioned QR SVD  ; this and 4
         are the most accurate options)
    4 :: CGEJSV (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) 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.
[in]LDX
    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 CGEEV.
    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 CGEEV for computing
    the eigenvalues of a Rayleigh quotient.
    If the call to CGEDMD 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_CGEEV),
    where LZWORK_CGEEV = MAX( 1, 2*N )  and the minimal
    LZWORK_SVD is calculated as follows
    If WHTSVD == 1 :: CGESVD ::
       LZWORK_SVD = MAX(1,2*MIN(M,N)+MAX(M,N))
    If WHTSVD == 2 :: CGESDD ::
       LZWORK_SVD = 2*MIN(M,N)*MIN(M,N)+2*MIN(M,N)+MAX(M,N)
    If WHTSVD == 3 :: CGESVDQ ::
       LZWORK_SVD = obtainable by a query
    If WHTSVD == 4 :: CGEJSV ::
       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 CGEDMD 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_CGEEV), where
    LRWORK_CGEEV = MAX(1,2*N) and RWORK_SVD is the real workspace
    for the SVD subroutine determined by the input parameter
    WHTSVD.
    If WHTSVD == 1 :: CGESVD ::
       LRWORK_SVD = 5*MIN(M,N)
    If WHTSVD == 2 :: CGESDD ::
       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 :: CGESVDQ ::
       LRWORK_SVD = obtainable by a query
    If WHTSVD == 4 :: CGEJSV ::
       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 496 of file cgedmd.f90.

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