Annotation of rpl/lapack/lapack/dsygv_2stage.f, revision 1.4

1.1       bertrand    1: *> \brief \b DSYGV_2STAGE
                      2: *
                      3: *  @precisions fortran d -> s
                      4: *
                      5: *  =========== DOCUMENTATION ===========
                      6: *
                      7: * Online html documentation available at
                      8: *            http://www.netlib.org/lapack/explore-html/
                      9: *
                     10: *> \htmlonly
                     11: *> Download DSYGV_2STAGE + dependencies
                     12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsygv_2stage.f">
                     13: *> [TGZ]</a>
                     14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsygv_2stage.f">
                     15: *> [ZIP]</a>
                     16: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsygv_2stage.f">
                     17: *> [TXT]</a>
                     18: *> \endhtmlonly
                     19: *
                     20: *  Definition:
                     21: *  ===========
                     22: *
                     23: *       SUBROUTINE DSYGV_2STAGE( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W,
                     24: *                                WORK, LWORK, INFO )
                     25: *
                     26: *       IMPLICIT NONE
                     27: *
                     28: *       .. Scalar Arguments ..
                     29: *       CHARACTER          JOBZ, UPLO
                     30: *       INTEGER            INFO, ITYPE, LDA, LDB, LWORK, N
                     31: *       ..
                     32: *       .. Array Arguments ..
                     33: *       DOUBLE PRECISION   A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
                     34: *       ..
                     35: *
                     36: *
                     37: *> \par Purpose:
                     38: *  =============
                     39: *>
                     40: *> \verbatim
                     41: *>
                     42: *> DSYGV_2STAGE computes all the eigenvalues, and optionally, the eigenvectors
                     43: *> of a real generalized symmetric-definite eigenproblem, of the form
                     44: *> A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.
                     45: *> Here A and B are assumed to be symmetric and B is also
                     46: *> positive definite.
                     47: *> This routine use the 2stage technique for the reduction to tridiagonal
                     48: *> which showed higher performance on recent architecture and for large
1.3       bertrand   49: *> sizes N>2000.
1.1       bertrand   50: *> \endverbatim
                     51: *
                     52: *  Arguments:
                     53: *  ==========
                     54: *
                     55: *> \param[in] ITYPE
                     56: *> \verbatim
                     57: *>          ITYPE is INTEGER
                     58: *>          Specifies the problem type to be solved:
                     59: *>          = 1:  A*x = (lambda)*B*x
                     60: *>          = 2:  A*B*x = (lambda)*x
                     61: *>          = 3:  B*A*x = (lambda)*x
                     62: *> \endverbatim
                     63: *>
                     64: *> \param[in] JOBZ
                     65: *> \verbatim
                     66: *>          JOBZ is CHARACTER*1
                     67: *>          = 'N':  Compute eigenvalues only;
                     68: *>          = 'V':  Compute eigenvalues and eigenvectors.
                     69: *>                  Not available in this release.
                     70: *> \endverbatim
                     71: *>
                     72: *> \param[in] UPLO
                     73: *> \verbatim
                     74: *>          UPLO is CHARACTER*1
                     75: *>          = 'U':  Upper triangles of A and B are stored;
                     76: *>          = 'L':  Lower triangles of A and B are stored.
                     77: *> \endverbatim
                     78: *>
                     79: *> \param[in] N
                     80: *> \verbatim
                     81: *>          N is INTEGER
                     82: *>          The order of the matrices A and B.  N >= 0.
                     83: *> \endverbatim
                     84: *>
                     85: *> \param[in,out] A
                     86: *> \verbatim
                     87: *>          A is DOUBLE PRECISION array, dimension (LDA, N)
                     88: *>          On entry, the symmetric matrix A.  If UPLO = 'U', the
                     89: *>          leading N-by-N upper triangular part of A contains the
                     90: *>          upper triangular part of the matrix A.  If UPLO = 'L',
                     91: *>          the leading N-by-N lower triangular part of A contains
                     92: *>          the lower triangular part of the matrix A.
                     93: *>
                     94: *>          On exit, if JOBZ = 'V', then if INFO = 0, A contains the
                     95: *>          matrix Z of eigenvectors.  The eigenvectors are normalized
                     96: *>          as follows:
                     97: *>          if ITYPE = 1 or 2, Z**T*B*Z = I;
                     98: *>          if ITYPE = 3, Z**T*inv(B)*Z = I.
                     99: *>          If JOBZ = 'N', then on exit the upper triangle (if UPLO='U')
                    100: *>          or the lower triangle (if UPLO='L') of A, including the
                    101: *>          diagonal, is destroyed.
                    102: *> \endverbatim
                    103: *>
                    104: *> \param[in] LDA
                    105: *> \verbatim
                    106: *>          LDA is INTEGER
                    107: *>          The leading dimension of the array A.  LDA >= max(1,N).
                    108: *> \endverbatim
                    109: *>
                    110: *> \param[in,out] B
                    111: *> \verbatim
                    112: *>          B is DOUBLE PRECISION array, dimension (LDB, N)
                    113: *>          On entry, the symmetric positive definite matrix B.
                    114: *>          If UPLO = 'U', the leading N-by-N upper triangular part of B
                    115: *>          contains the upper triangular part of the matrix B.
                    116: *>          If UPLO = 'L', the leading N-by-N lower triangular part of B
                    117: *>          contains the lower triangular part of the matrix B.
                    118: *>
                    119: *>          On exit, if INFO <= N, the part of B containing the matrix is
                    120: *>          overwritten by the triangular factor U or L from the Cholesky
                    121: *>          factorization B = U**T*U or B = L*L**T.
                    122: *> \endverbatim
                    123: *>
                    124: *> \param[in] LDB
                    125: *> \verbatim
                    126: *>          LDB is INTEGER
                    127: *>          The leading dimension of the array B.  LDB >= max(1,N).
                    128: *> \endverbatim
                    129: *>
                    130: *> \param[out] W
                    131: *> \verbatim
                    132: *>          W is DOUBLE PRECISION array, dimension (N)
                    133: *>          If INFO = 0, the eigenvalues in ascending order.
                    134: *> \endverbatim
                    135: *>
                    136: *> \param[out] WORK
                    137: *> \verbatim
                    138: *>          WORK is DOUBLE PRECISION array, dimension (MAX(1,LWORK))
                    139: *>          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
                    140: *> \endverbatim
                    141: *>
                    142: *> \param[in] LWORK
                    143: *> \verbatim
                    144: *>          LWORK is INTEGER
                    145: *>          The length of the array WORK. LWORK >= 1, when N <= 1;
                    146: *>          otherwise  
                    147: *>          If JOBZ = 'N' and N > 1, LWORK must be queried.
                    148: *>                                   LWORK = MAX(1, dimension) where
                    149: *>                                   dimension = max(stage1,stage2) + (KD+1)*N + 2*N
                    150: *>                                             = N*KD + N*max(KD+1,FACTOPTNB) 
                    151: *>                                               + max(2*KD*KD, KD*NTHREADS) 
                    152: *>                                               + (KD+1)*N + 2*N
                    153: *>                                   where KD is the blocking size of the reduction,
                    154: *>                                   FACTOPTNB is the blocking used by the QR or LQ
                    155: *>                                   algorithm, usually FACTOPTNB=128 is a good choice
                    156: *>                                   NTHREADS is the number of threads used when
                    157: *>                                   openMP compilation is enabled, otherwise =1.
                    158: *>          If JOBZ = 'V' and N > 1, LWORK must be queried. Not yet available
                    159: *>
                    160: *>          If LWORK = -1, then a workspace query is assumed; the routine
                    161: *>          only calculates the optimal size of the WORK array, returns
                    162: *>          this value as the first entry of the WORK array, and no error
                    163: *>          message related to LWORK is issued by XERBLA.
                    164: *> \endverbatim
                    165: *>
                    166: *> \param[out] INFO
                    167: *> \verbatim
                    168: *>          INFO is INTEGER
                    169: *>          = 0:  successful exit
                    170: *>          < 0:  if INFO = -i, the i-th argument had an illegal value
                    171: *>          > 0:  DPOTRF or DSYEV returned an error code:
                    172: *>             <= N:  if INFO = i, DSYEV failed to converge;
                    173: *>                    i off-diagonal elements of an intermediate
                    174: *>                    tridiagonal form did not converge to zero;
                    175: *>             > N:   if INFO = N + i, for 1 <= i <= N, then the leading
                    176: *>                    minor of order i of B is not positive definite.
                    177: *>                    The factorization of B could not be completed and
                    178: *>                    no eigenvalues or eigenvectors were computed.
                    179: *> \endverbatim
                    180: *
                    181: *  Authors:
                    182: *  ========
                    183: *
                    184: *> \author Univ. of Tennessee
                    185: *> \author Univ. of California Berkeley
                    186: *> \author Univ. of Colorado Denver
                    187: *> \author NAG Ltd.
                    188: *
1.3       bertrand  189: *> \date November 2017
1.1       bertrand  190: *
                    191: *> \ingroup doubleSYeigen
                    192: *
                    193: *> \par Further Details:
                    194: *  =====================
                    195: *>
                    196: *> \verbatim
                    197: *>
                    198: *>  All details about the 2stage techniques are available in:
                    199: *>
                    200: *>  Azzam Haidar, Hatem Ltaief, and Jack Dongarra.
                    201: *>  Parallel reduction to condensed forms for symmetric eigenvalue problems
                    202: *>  using aggregated fine-grained and memory-aware kernels. In Proceedings
                    203: *>  of 2011 International Conference for High Performance Computing,
                    204: *>  Networking, Storage and Analysis (SC '11), New York, NY, USA,
                    205: *>  Article 8 , 11 pages.
                    206: *>  http://doi.acm.org/10.1145/2063384.2063394
                    207: *>
                    208: *>  A. Haidar, J. Kurzak, P. Luszczek, 2013.
                    209: *>  An improved parallel singular value algorithm and its implementation 
                    210: *>  for multicore hardware, In Proceedings of 2013 International Conference
                    211: *>  for High Performance Computing, Networking, Storage and Analysis (SC '13).
                    212: *>  Denver, Colorado, USA, 2013.
                    213: *>  Article 90, 12 pages.
                    214: *>  http://doi.acm.org/10.1145/2503210.2503292
                    215: *>
                    216: *>  A. Haidar, R. Solca, S. Tomov, T. Schulthess and J. Dongarra.
                    217: *>  A novel hybrid CPU-GPU generalized eigensolver for electronic structure 
                    218: *>  calculations based on fine-grained memory aware tasks.
                    219: *>  International Journal of High Performance Computing Applications.
                    220: *>  Volume 28 Issue 2, Pages 196-209, May 2014.
                    221: *>  http://hpc.sagepub.com/content/28/2/196 
                    222: *>
                    223: *> \endverbatim
                    224: *
                    225: *  =====================================================================
                    226:       SUBROUTINE DSYGV_2STAGE( ITYPE, JOBZ, UPLO, N, A, LDA, B, LDB, W,
                    227:      $                         WORK, LWORK, INFO )
                    228: *
                    229:       IMPLICIT NONE
                    230: *
1.3       bertrand  231: *  -- LAPACK driver routine (version 3.8.0) --
1.1       bertrand  232: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
                    233: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
1.3       bertrand  234: *     November 2017
1.1       bertrand  235: *
                    236: *     .. Scalar Arguments ..
                    237:       CHARACTER          JOBZ, UPLO
                    238:       INTEGER            INFO, ITYPE, LDA, LDB, LWORK, N
                    239: *     ..
                    240: *     .. Array Arguments ..
                    241:       DOUBLE PRECISION   A( LDA, * ), B( LDB, * ), W( * ), WORK( * )
                    242: *     ..
                    243: *
                    244: *  =====================================================================
                    245: *
                    246: *     .. Parameters ..
                    247:       DOUBLE PRECISION   ONE
                    248:       PARAMETER          ( ONE = 1.0D+0 )
                    249: *     ..
                    250: *     .. Local Scalars ..
                    251:       LOGICAL            LQUERY, UPPER, WANTZ
                    252:       CHARACTER          TRANS
                    253:       INTEGER            NEIG, LWMIN, LHTRD, LWTRD, KD, IB 
                    254: *     ..
                    255: *     .. External Functions ..
                    256:       LOGICAL            LSAME
1.3       bertrand  257:       INTEGER            ILAENV2STAGE
                    258:       EXTERNAL           LSAME, ILAENV2STAGE
1.1       bertrand  259: *     ..
                    260: *     .. External Subroutines ..
                    261:       EXTERNAL           DPOTRF, DSYGST, DTRMM, DTRSM, XERBLA,
                    262:      $                   DSYEV_2STAGE
                    263: *     ..
                    264: *     .. Intrinsic Functions ..
                    265:       INTRINSIC          MAX
                    266: *     ..
                    267: *     .. Executable Statements ..
                    268: *
                    269: *     Test the input parameters.
                    270: *
                    271:       WANTZ = LSAME( JOBZ, 'V' )
                    272:       UPPER = LSAME( UPLO, 'U' )
                    273:       LQUERY = ( LWORK.EQ.-1 )
                    274: *
                    275:       INFO = 0
                    276:       IF( ITYPE.LT.1 .OR. ITYPE.GT.3 ) THEN
                    277:          INFO = -1
                    278:       ELSE IF( .NOT.( LSAME( JOBZ, 'N' ) ) ) THEN
                    279:          INFO = -2
                    280:       ELSE IF( .NOT.( UPPER .OR. LSAME( UPLO, 'L' ) ) ) THEN
                    281:          INFO = -3
                    282:       ELSE IF( N.LT.0 ) THEN
                    283:          INFO = -4
                    284:       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
                    285:          INFO = -6
                    286:       ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
                    287:          INFO = -8
                    288:       END IF
                    289: *
                    290:       IF( INFO.EQ.0 ) THEN
1.3       bertrand  291:          KD    = ILAENV2STAGE( 1, 'DSYTRD_2STAGE', JOBZ, N, -1, -1, -1 )
                    292:          IB    = ILAENV2STAGE( 2, 'DSYTRD_2STAGE', JOBZ, N, KD, -1, -1 )
                    293:          LHTRD = ILAENV2STAGE( 3, 'DSYTRD_2STAGE', JOBZ, N, KD, IB, -1 )
                    294:          LWTRD = ILAENV2STAGE( 4, 'DSYTRD_2STAGE', JOBZ, N, KD, IB, -1 )
1.1       bertrand  295:          LWMIN = 2*N + LHTRD + LWTRD
                    296:          WORK( 1 )  = LWMIN
                    297: *
                    298:          IF( LWORK.LT.LWMIN .AND. .NOT.LQUERY ) THEN
                    299:             INFO = -11
                    300:          END IF
                    301:       END IF
                    302: *
                    303:       IF( INFO.NE.0 ) THEN
                    304:          CALL XERBLA( 'DSYGV_2STAGE ', -INFO )
                    305:          RETURN
                    306:       ELSE IF( LQUERY ) THEN
                    307:          RETURN
                    308:       END IF
                    309: *
                    310: *     Quick return if possible
                    311: *
                    312:       IF( N.EQ.0 )
                    313:      $   RETURN
                    314: *
                    315: *     Form a Cholesky factorization of B.
                    316: *
                    317:       CALL DPOTRF( UPLO, N, B, LDB, INFO )
                    318:       IF( INFO.NE.0 ) THEN
                    319:          INFO = N + INFO
                    320:          RETURN
                    321:       END IF
                    322: *
                    323: *     Transform problem to standard eigenvalue problem and solve.
                    324: *
                    325:       CALL DSYGST( ITYPE, UPLO, N, A, LDA, B, LDB, INFO )
                    326:       CALL DSYEV_2STAGE( JOBZ, UPLO, N, A, LDA, W, WORK, LWORK, INFO )
                    327: *
                    328:       IF( WANTZ ) THEN
                    329: *
                    330: *        Backtransform eigenvectors to the original problem.
                    331: *
                    332:          NEIG = N
                    333:          IF( INFO.GT.0 )
                    334:      $      NEIG = INFO - 1
                    335:          IF( ITYPE.EQ.1 .OR. ITYPE.EQ.2 ) THEN
                    336: *
                    337: *           For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
                    338: *           backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
                    339: *
                    340:             IF( UPPER ) THEN
                    341:                TRANS = 'N'
                    342:             ELSE
                    343:                TRANS = 'T'
                    344:             END IF
                    345: *
                    346:             CALL DTRSM( 'Left', UPLO, TRANS, 'Non-unit', N, NEIG, ONE,
                    347:      $                  B, LDB, A, LDA )
                    348: *
                    349:          ELSE IF( ITYPE.EQ.3 ) THEN
                    350: *
                    351: *           For B*A*x=(lambda)*x;
                    352: *           backtransform eigenvectors: x = L*y or U**T*y
                    353: *
                    354:             IF( UPPER ) THEN
                    355:                TRANS = 'T'
                    356:             ELSE
                    357:                TRANS = 'N'
                    358:             END IF
                    359: *
                    360:             CALL DTRMM( 'Left', UPLO, TRANS, 'Non-unit', N, NEIG, ONE,
                    361:      $                  B, LDB, A, LDA )
                    362:          END IF
                    363:       END IF
                    364: *
                    365:       WORK( 1 ) = LWMIN
                    366:       RETURN
                    367: *
                    368: *     End of DSYGV_2STAGE
                    369: *
                    370:       END

CVSweb interface <joel.bertrand@systella.fr>