Annotation of rpl/lapack/lapack/dsysvxx.f, revision 1.5

1.5     ! bertrand    1: *> \brief \b DSYSVXX
        !             2: *
        !             3: *  =========== DOCUMENTATION ===========
        !             4: *
        !             5: * Online html documentation available at 
        !             6: *            http://www.netlib.org/lapack/explore-html/ 
        !             7: *
        !             8: *> \htmlonly
        !             9: *> Download DSYSVXX + dependencies 
        !            10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dsysvxx.f"> 
        !            11: *> [TGZ]</a> 
        !            12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dsysvxx.f"> 
        !            13: *> [ZIP]</a> 
        !            14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dsysvxx.f"> 
        !            15: *> [TXT]</a>
        !            16: *> \endhtmlonly 
        !            17: *
        !            18: *  Definition:
        !            19: *  ===========
        !            20: *
        !            21: *       SUBROUTINE DSYSVXX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, IPIV,
        !            22: *                           EQUED, S, B, LDB, X, LDX, RCOND, RPVGRW, BERR,
        !            23: *                           N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
        !            24: *                           NPARAMS, PARAMS, WORK, IWORK, INFO )
        !            25: * 
        !            26: *       .. Scalar Arguments ..
        !            27: *       CHARACTER          EQUED, FACT, UPLO
        !            28: *       INTEGER            INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
        !            29: *      $                   N_ERR_BNDS
        !            30: *       DOUBLE PRECISION   RCOND, RPVGRW
        !            31: *       ..
        !            32: *       .. Array Arguments ..
        !            33: *       INTEGER            IPIV( * ), IWORK( * )
        !            34: *       DOUBLE PRECISION   A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
        !            35: *      $                   X( LDX, * ), WORK( * )
        !            36: *       DOUBLE PRECISION   S( * ), PARAMS( * ), BERR( * ),
        !            37: *      $                   ERR_BNDS_NORM( NRHS, * ),
        !            38: *      $                   ERR_BNDS_COMP( NRHS, * )
        !            39: *       ..
        !            40: *  
        !            41: *
        !            42: *> \par Purpose:
        !            43: *  =============
        !            44: *>
        !            45: *> \verbatim
        !            46: *>
        !            47: *>    DSYSVXX uses the diagonal pivoting factorization to compute the
        !            48: *>    solution to a double precision system of linear equations A * X = B, where A
        !            49: *>    is an N-by-N symmetric matrix and X and B are N-by-NRHS matrices.
        !            50: *>
        !            51: *>    If requested, both normwise and maximum componentwise error bounds
        !            52: *>    are returned. DSYSVXX will return a solution with a tiny
        !            53: *>    guaranteed error (O(eps) where eps is the working machine
        !            54: *>    precision) unless the matrix is very ill-conditioned, in which
        !            55: *>    case a warning is returned. Relevant condition numbers also are
        !            56: *>    calculated and returned.
        !            57: *>
        !            58: *>    DSYSVXX accepts user-provided factorizations and equilibration
        !            59: *>    factors; see the definitions of the FACT and EQUED options.
        !            60: *>    Solving with refinement and using a factorization from a previous
        !            61: *>    DSYSVXX call will also produce a solution with either O(eps)
        !            62: *>    errors or warnings, but we cannot make that claim for general
        !            63: *>    user-provided factorizations and equilibration factors if they
        !            64: *>    differ from what DSYSVXX would itself produce.
        !            65: *> \endverbatim
        !            66: *
        !            67: *> \par Description:
        !            68: *  =================
        !            69: *>
        !            70: *> \verbatim
        !            71: *>
        !            72: *>    The following steps are performed:
        !            73: *>
        !            74: *>    1. If FACT = 'E', double precision scaling factors are computed to equilibrate
        !            75: *>    the system:
        !            76: *>
        !            77: *>      diag(S)*A*diag(S)     *inv(diag(S))*X = diag(S)*B
        !            78: *>
        !            79: *>    Whether or not the system will be equilibrated depends on the
        !            80: *>    scaling of the matrix A, but if equilibration is used, A is
        !            81: *>    overwritten by diag(S)*A*diag(S) and B by diag(S)*B.
        !            82: *>
        !            83: *>    2. If FACT = 'N' or 'E', the LU decomposition is used to factor
        !            84: *>    the matrix A (after equilibration if FACT = 'E') as
        !            85: *>
        !            86: *>       A = U * D * U**T,  if UPLO = 'U', or
        !            87: *>       A = L * D * L**T,  if UPLO = 'L',
        !            88: *>
        !            89: *>    where U (or L) is a product of permutation and unit upper (lower)
        !            90: *>    triangular matrices, and D is symmetric and block diagonal with
        !            91: *>    1-by-1 and 2-by-2 diagonal blocks.
        !            92: *>
        !            93: *>    3. If some D(i,i)=0, so that D is exactly singular, then the
        !            94: *>    routine returns with INFO = i. Otherwise, the factored form of A
        !            95: *>    is used to estimate the condition number of the matrix A (see
        !            96: *>    argument RCOND).  If the reciprocal of the condition number is
        !            97: *>    less than machine precision, the routine still goes on to solve
        !            98: *>    for X and compute error bounds as described below.
        !            99: *>
        !           100: *>    4. The system of equations is solved for X using the factored form
        !           101: *>    of A.
        !           102: *>
        !           103: *>    5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
        !           104: *>    the routine will use iterative refinement to try to get a small
        !           105: *>    error and error bounds.  Refinement calculates the residual to at
        !           106: *>    least twice the working precision.
        !           107: *>
        !           108: *>    6. If equilibration was used, the matrix X is premultiplied by
        !           109: *>    diag(R) so that it solves the original system before
        !           110: *>    equilibration.
        !           111: *> \endverbatim
        !           112: *
        !           113: *  Arguments:
        !           114: *  ==========
        !           115: *
        !           116: *> \verbatim
        !           117: *>     Some optional parameters are bundled in the PARAMS array.  These
        !           118: *>     settings determine how refinement is performed, but often the
        !           119: *>     defaults are acceptable.  If the defaults are acceptable, users
        !           120: *>     can pass NPARAMS = 0 which prevents the source code from accessing
        !           121: *>     the PARAMS argument.
        !           122: *> \endverbatim
        !           123: *>
        !           124: *> \param[in] FACT
        !           125: *> \verbatim
        !           126: *>          FACT is CHARACTER*1
        !           127: *>     Specifies whether or not the factored form of the matrix A is
        !           128: *>     supplied on entry, and if not, whether the matrix A should be
        !           129: *>     equilibrated before it is factored.
        !           130: *>       = 'F':  On entry, AF and IPIV contain the factored form of A.
        !           131: *>               If EQUED is not 'N', the matrix A has been
        !           132: *>               equilibrated with scaling factors given by S.
        !           133: *>               A, AF, and IPIV are not modified.
        !           134: *>       = 'N':  The matrix A will be copied to AF and factored.
        !           135: *>       = 'E':  The matrix A will be equilibrated if necessary, then
        !           136: *>               copied to AF and factored.
        !           137: *> \endverbatim
        !           138: *>
        !           139: *> \param[in] UPLO
        !           140: *> \verbatim
        !           141: *>          UPLO is CHARACTER*1
        !           142: *>       = 'U':  Upper triangle of A is stored;
        !           143: *>       = 'L':  Lower triangle of A is stored.
        !           144: *> \endverbatim
        !           145: *>
        !           146: *> \param[in] N
        !           147: *> \verbatim
        !           148: *>          N is INTEGER
        !           149: *>     The number of linear equations, i.e., the order of the
        !           150: *>     matrix A.  N >= 0.
        !           151: *> \endverbatim
        !           152: *>
        !           153: *> \param[in] NRHS
        !           154: *> \verbatim
        !           155: *>          NRHS is INTEGER
        !           156: *>     The number of right hand sides, i.e., the number of columns
        !           157: *>     of the matrices B and X.  NRHS >= 0.
        !           158: *> \endverbatim
        !           159: *>
        !           160: *> \param[in,out] A
        !           161: *> \verbatim
        !           162: *>          A is DOUBLE PRECISION array, dimension (LDA,N)
        !           163: *>     The symmetric matrix A.  If UPLO = 'U', the leading N-by-N
        !           164: *>     upper triangular part of A contains the upper triangular
        !           165: *>     part of the matrix A, and the strictly lower triangular
        !           166: *>     part of A is not referenced.  If UPLO = 'L', the leading
        !           167: *>     N-by-N lower triangular part of A contains the lower
        !           168: *>     triangular part of the matrix A, and the strictly upper
        !           169: *>     triangular part of A is not referenced.
        !           170: *>
        !           171: *>     On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
        !           172: *>     diag(S)*A*diag(S).
        !           173: *> \endverbatim
        !           174: *>
        !           175: *> \param[in] LDA
        !           176: *> \verbatim
        !           177: *>          LDA is INTEGER
        !           178: *>     The leading dimension of the array A.  LDA >= max(1,N).
        !           179: *> \endverbatim
        !           180: *>
        !           181: *> \param[in,out] AF
        !           182: *> \verbatim
        !           183: *>          AF is or output) DOUBLE PRECISION array, dimension (LDAF,N)
        !           184: *>     If FACT = 'F', then AF is an input argument and on entry
        !           185: *>     contains the block diagonal matrix D and the multipliers
        !           186: *>     used to obtain the factor U or L from the factorization A =
        !           187: *>     U*D*U**T or A = L*D*L**T as computed by DSYTRF.
        !           188: *>
        !           189: *>     If FACT = 'N', then AF is an output argument and on exit
        !           190: *>     returns the block diagonal matrix D and the multipliers
        !           191: *>     used to obtain the factor U or L from the factorization A =
        !           192: *>     U*D*U**T or A = L*D*L**T.
        !           193: *> \endverbatim
        !           194: *>
        !           195: *> \param[in] LDAF
        !           196: *> \verbatim
        !           197: *>          LDAF is INTEGER
        !           198: *>     The leading dimension of the array AF.  LDAF >= max(1,N).
        !           199: *> \endverbatim
        !           200: *>
        !           201: *> \param[in,out] IPIV
        !           202: *> \verbatim
        !           203: *>          IPIV is or output) INTEGER array, dimension (N)
        !           204: *>     If FACT = 'F', then IPIV is an input argument and on entry
        !           205: *>     contains details of the interchanges and the block
        !           206: *>     structure of D, as determined by DSYTRF.  If IPIV(k) > 0,
        !           207: *>     then rows and columns k and IPIV(k) were interchanged and
        !           208: *>     D(k,k) is a 1-by-1 diagonal block.  If UPLO = 'U' and
        !           209: *>     IPIV(k) = IPIV(k-1) < 0, then rows and columns k-1 and
        !           210: *>     -IPIV(k) were interchanged and D(k-1:k,k-1:k) is a 2-by-2
        !           211: *>     diagonal block.  If UPLO = 'L' and IPIV(k) = IPIV(k+1) < 0,
        !           212: *>     then rows and columns k+1 and -IPIV(k) were interchanged
        !           213: *>     and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
        !           214: *>
        !           215: *>     If FACT = 'N', then IPIV is an output argument and on exit
        !           216: *>     contains details of the interchanges and the block
        !           217: *>     structure of D, as determined by DSYTRF.
        !           218: *> \endverbatim
        !           219: *>
        !           220: *> \param[in,out] EQUED
        !           221: *> \verbatim
        !           222: *>          EQUED is or output) CHARACTER*1
        !           223: *>     Specifies the form of equilibration that was done.
        !           224: *>       = 'N':  No equilibration (always true if FACT = 'N').
        !           225: *>       = 'Y':  Both row and column equilibration, i.e., A has been
        !           226: *>               replaced by diag(S) * A * diag(S).
        !           227: *>     EQUED is an input argument if FACT = 'F'; otherwise, it is an
        !           228: *>     output argument.
        !           229: *> \endverbatim
        !           230: *>
        !           231: *> \param[in,out] S
        !           232: *> \verbatim
        !           233: *>          S is or output) DOUBLE PRECISION array, dimension (N)
        !           234: *>     The scale factors for A.  If EQUED = 'Y', A is multiplied on
        !           235: *>     the left and right by diag(S).  S is an input argument if FACT =
        !           236: *>     'F'; otherwise, S is an output argument.  If FACT = 'F' and EQUED
        !           237: *>     = 'Y', each element of S must be positive.  If S is output, each
        !           238: *>     element of S is a power of the radix. If S is input, each element
        !           239: *>     of S should be a power of the radix to ensure a reliable solution
        !           240: *>     and error estimates. Scaling by powers of the radix does not cause
        !           241: *>     rounding errors unless the result underflows or overflows.
        !           242: *>     Rounding errors during scaling lead to refining with a matrix that
        !           243: *>     is not equivalent to the input matrix, producing error estimates
        !           244: *>     that may not be reliable.
        !           245: *> \endverbatim
        !           246: *>
        !           247: *> \param[in,out] B
        !           248: *> \verbatim
        !           249: *>          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
        !           250: *>     On entry, the N-by-NRHS right hand side matrix B.
        !           251: *>     On exit,
        !           252: *>     if EQUED = 'N', B is not modified;
        !           253: *>     if EQUED = 'Y', B is overwritten by diag(S)*B;
        !           254: *> \endverbatim
        !           255: *>
        !           256: *> \param[in] LDB
        !           257: *> \verbatim
        !           258: *>          LDB is INTEGER
        !           259: *>     The leading dimension of the array B.  LDB >= max(1,N).
        !           260: *> \endverbatim
        !           261: *>
        !           262: *> \param[out] X
        !           263: *> \verbatim
        !           264: *>          X is DOUBLE PRECISION array, dimension (LDX,NRHS)
        !           265: *>     If INFO = 0, the N-by-NRHS solution matrix X to the original
        !           266: *>     system of equations.  Note that A and B are modified on exit if
        !           267: *>     EQUED .ne. 'N', and the solution to the equilibrated system is
        !           268: *>     inv(diag(S))*X.
        !           269: *> \endverbatim
        !           270: *>
        !           271: *> \param[in] LDX
        !           272: *> \verbatim
        !           273: *>          LDX is INTEGER
        !           274: *>     The leading dimension of the array X.  LDX >= max(1,N).
        !           275: *> \endverbatim
        !           276: *>
        !           277: *> \param[out] RCOND
        !           278: *> \verbatim
        !           279: *>          RCOND is DOUBLE PRECISION
        !           280: *>     Reciprocal scaled condition number.  This is an estimate of the
        !           281: *>     reciprocal Skeel condition number of the matrix A after
        !           282: *>     equilibration (if done).  If this is less than the machine
        !           283: *>     precision (in particular, if it is zero), the matrix is singular
        !           284: *>     to working precision.  Note that the error may still be small even
        !           285: *>     if this number is very small and the matrix appears ill-
        !           286: *>     conditioned.
        !           287: *> \endverbatim
        !           288: *>
        !           289: *> \param[out] RPVGRW
        !           290: *> \verbatim
        !           291: *>          RPVGRW is DOUBLE PRECISION
        !           292: *>     Reciprocal pivot growth.  On exit, this contains the reciprocal
        !           293: *>     pivot growth factor norm(A)/norm(U). The "max absolute element"
        !           294: *>     norm is used.  If this is much less than 1, then the stability of
        !           295: *>     the LU factorization of the (equilibrated) matrix A could be poor.
        !           296: *>     This also means that the solution X, estimated condition numbers,
        !           297: *>     and error bounds could be unreliable. If factorization fails with
        !           298: *>     0<INFO<=N, then this contains the reciprocal pivot growth factor
        !           299: *>     for the leading INFO columns of A.
        !           300: *> \endverbatim
        !           301: *>
        !           302: *> \param[out] BERR
        !           303: *> \verbatim
        !           304: *>          BERR is DOUBLE PRECISION array, dimension (NRHS)
        !           305: *>     Componentwise relative backward error.  This is the
        !           306: *>     componentwise relative backward error of each solution vector X(j)
        !           307: *>     (i.e., the smallest relative change in any element of A or B that
        !           308: *>     makes X(j) an exact solution).
        !           309: *> \endverbatim
        !           310: *>
        !           311: *> \param[in] N_ERR_BNDS
        !           312: *> \verbatim
        !           313: *>          N_ERR_BNDS is INTEGER
        !           314: *>     Number of error bounds to return for each right hand side
        !           315: *>     and each type (normwise or componentwise).  See ERR_BNDS_NORM and
        !           316: *>     ERR_BNDS_COMP below.
        !           317: *> \endverbatim
        !           318: *>
        !           319: *> \param[out] ERR_BNDS_NORM
        !           320: *> \verbatim
        !           321: *>          ERR_BNDS_NORM is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
        !           322: *>     For each right-hand side, this array contains information about
        !           323: *>     various error bounds and condition numbers corresponding to the
        !           324: *>     normwise relative error, which is defined as follows:
        !           325: *>
        !           326: *>     Normwise relative error in the ith solution vector:
        !           327: *>             max_j (abs(XTRUE(j,i) - X(j,i)))
        !           328: *>            ------------------------------
        !           329: *>                  max_j abs(X(j,i))
        !           330: *>
        !           331: *>     The array is indexed by the type of error information as described
        !           332: *>     below. There currently are up to three pieces of information
        !           333: *>     returned.
        !           334: *>
        !           335: *>     The first index in ERR_BNDS_NORM(i,:) corresponds to the ith
        !           336: *>     right-hand side.
        !           337: *>
        !           338: *>     The second index in ERR_BNDS_NORM(:,err) contains the following
        !           339: *>     three fields:
        !           340: *>     err = 1 "Trust/don't trust" boolean. Trust the answer if the
        !           341: *>              reciprocal condition number is less than the threshold
        !           342: *>              sqrt(n) * dlamch('Epsilon').
        !           343: *>
        !           344: *>     err = 2 "Guaranteed" error bound: The estimated forward error,
        !           345: *>              almost certainly within a factor of 10 of the true error
        !           346: *>              so long as the next entry is greater than the threshold
        !           347: *>              sqrt(n) * dlamch('Epsilon'). This error bound should only
        !           348: *>              be trusted if the previous boolean is true.
        !           349: *>
        !           350: *>     err = 3  Reciprocal condition number: Estimated normwise
        !           351: *>              reciprocal condition number.  Compared with the threshold
        !           352: *>              sqrt(n) * dlamch('Epsilon') to determine if the error
        !           353: *>              estimate is "guaranteed". These reciprocal condition
        !           354: *>              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
        !           355: *>              appropriately scaled matrix Z.
        !           356: *>              Let Z = S*A, where S scales each row by a power of the
        !           357: *>              radix so all absolute row sums of Z are approximately 1.
        !           358: *>
        !           359: *>     See Lapack Working Note 165 for further details and extra
        !           360: *>     cautions.
        !           361: *> \endverbatim
        !           362: *>
        !           363: *> \param[out] ERR_BNDS_COMP
        !           364: *> \verbatim
        !           365: *>          ERR_BNDS_COMP is DOUBLE PRECISION array, dimension (NRHS, N_ERR_BNDS)
        !           366: *>     For each right-hand side, this array contains information about
        !           367: *>     various error bounds and condition numbers corresponding to the
        !           368: *>     componentwise relative error, which is defined as follows:
        !           369: *>
        !           370: *>     Componentwise relative error in the ith solution vector:
        !           371: *>                    abs(XTRUE(j,i) - X(j,i))
        !           372: *>             max_j ----------------------
        !           373: *>                         abs(X(j,i))
        !           374: *>
        !           375: *>     The array is indexed by the right-hand side i (on which the
        !           376: *>     componentwise relative error depends), and the type of error
        !           377: *>     information as described below. There currently are up to three
        !           378: *>     pieces of information returned for each right-hand side. If
        !           379: *>     componentwise accuracy is not requested (PARAMS(3) = 0.0), then
        !           380: *>     ERR_BNDS_COMP is not accessed.  If N_ERR_BNDS .LT. 3, then at most
        !           381: *>     the first (:,N_ERR_BNDS) entries are returned.
        !           382: *>
        !           383: *>     The first index in ERR_BNDS_COMP(i,:) corresponds to the ith
        !           384: *>     right-hand side.
        !           385: *>
        !           386: *>     The second index in ERR_BNDS_COMP(:,err) contains the following
        !           387: *>     three fields:
        !           388: *>     err = 1 "Trust/don't trust" boolean. Trust the answer if the
        !           389: *>              reciprocal condition number is less than the threshold
        !           390: *>              sqrt(n) * dlamch('Epsilon').
        !           391: *>
        !           392: *>     err = 2 "Guaranteed" error bound: The estimated forward error,
        !           393: *>              almost certainly within a factor of 10 of the true error
        !           394: *>              so long as the next entry is greater than the threshold
        !           395: *>              sqrt(n) * dlamch('Epsilon'). This error bound should only
        !           396: *>              be trusted if the previous boolean is true.
        !           397: *>
        !           398: *>     err = 3  Reciprocal condition number: Estimated componentwise
        !           399: *>              reciprocal condition number.  Compared with the threshold
        !           400: *>              sqrt(n) * dlamch('Epsilon') to determine if the error
        !           401: *>              estimate is "guaranteed". These reciprocal condition
        !           402: *>              numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for some
        !           403: *>              appropriately scaled matrix Z.
        !           404: *>              Let Z = S*(A*diag(x)), where x is the solution for the
        !           405: *>              current right-hand side and S scales each row of
        !           406: *>              A*diag(x) by a power of the radix so all absolute row
        !           407: *>              sums of Z are approximately 1.
        !           408: *>
        !           409: *>     See Lapack Working Note 165 for further details and extra
        !           410: *>     cautions.
        !           411: *> \endverbatim
        !           412: *>
        !           413: *> \param[in] NPARAMS
        !           414: *> \verbatim
        !           415: *>          NPARAMS is INTEGER
        !           416: *>     Specifies the number of parameters set in PARAMS.  If .LE. 0, the
        !           417: *>     PARAMS array is never referenced and default values are used.
        !           418: *> \endverbatim
        !           419: *>
        !           420: *> \param[in,out] PARAMS
        !           421: *> \verbatim
        !           422: *>          PARAMS is / output) DOUBLE PRECISION array, dimension (NPARAMS)
        !           423: *>     Specifies algorithm parameters.  If an entry is .LT. 0.0, then
        !           424: *>     that entry will be filled with default value used for that
        !           425: *>     parameter.  Only positions up to NPARAMS are accessed; defaults
        !           426: *>     are used for higher-numbered parameters.
        !           427: *>
        !           428: *>       PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
        !           429: *>            refinement or not.
        !           430: *>         Default: 1.0D+0
        !           431: *>            = 0.0 : No refinement is performed, and no error bounds are
        !           432: *>                    computed.
        !           433: *>            = 1.0 : Use the extra-precise refinement algorithm.
        !           434: *>              (other values are reserved for future use)
        !           435: *>
        !           436: *>       PARAMS(LA_LINRX_ITHRESH_I = 2) : Maximum number of residual
        !           437: *>            computations allowed for refinement.
        !           438: *>         Default: 10
        !           439: *>         Aggressive: Set to 100 to permit convergence using approximate
        !           440: *>                     factorizations or factorizations other than LU. If
        !           441: *>                     the factorization uses a technique other than
        !           442: *>                     Gaussian elimination, the guarantees in
        !           443: *>                     err_bnds_norm and err_bnds_comp may no longer be
        !           444: *>                     trustworthy.
        !           445: *>
        !           446: *>       PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining if the code
        !           447: *>            will attempt to find a solution with small componentwise
        !           448: *>            relative error in the double-precision algorithm.  Positive
        !           449: *>            is true, 0.0 is false.
        !           450: *>         Default: 1.0 (attempt componentwise convergence)
        !           451: *> \endverbatim
        !           452: *>
        !           453: *> \param[out] WORK
        !           454: *> \verbatim
        !           455: *>          WORK is DOUBLE PRECISION array, dimension (4*N)
        !           456: *> \endverbatim
        !           457: *>
        !           458: *> \param[out] IWORK
        !           459: *> \verbatim
        !           460: *>          IWORK is INTEGER array, dimension (N)
        !           461: *> \endverbatim
        !           462: *>
        !           463: *> \param[out] INFO
        !           464: *> \verbatim
        !           465: *>          INFO is INTEGER
        !           466: *>       = 0:  Successful exit. The solution to every right-hand side is
        !           467: *>         guaranteed.
        !           468: *>       < 0:  If INFO = -i, the i-th argument had an illegal value
        !           469: *>       > 0 and <= N:  U(INFO,INFO) is exactly zero.  The factorization
        !           470: *>         has been completed, but the factor U is exactly singular, so
        !           471: *>         the solution and error bounds could not be computed. RCOND = 0
        !           472: *>         is returned.
        !           473: *>       = N+J: The solution corresponding to the Jth right-hand side is
        !           474: *>         not guaranteed. The solutions corresponding to other right-
        !           475: *>         hand sides K with K > J may not be guaranteed as well, but
        !           476: *>         only the first such right-hand side is reported. If a small
        !           477: *>         componentwise error is not requested (PARAMS(3) = 0.0) then
        !           478: *>         the Jth right-hand side is the first with a normwise error
        !           479: *>         bound that is not guaranteed (the smallest J such
        !           480: *>         that ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0)
        !           481: *>         the Jth right-hand side is the first with either a normwise or
        !           482: *>         componentwise error bound that is not guaranteed (the smallest
        !           483: *>         J such that either ERR_BNDS_NORM(J,1) = 0.0 or
        !           484: *>         ERR_BNDS_COMP(J,1) = 0.0). See the definition of
        !           485: *>         ERR_BNDS_NORM(:,1) and ERR_BNDS_COMP(:,1). To get information
        !           486: *>         about all of the right-hand sides check ERR_BNDS_NORM or
        !           487: *>         ERR_BNDS_COMP.
        !           488: *> \endverbatim
        !           489: *
        !           490: *  Authors:
        !           491: *  ========
        !           492: *
        !           493: *> \author Univ. of Tennessee 
        !           494: *> \author Univ. of California Berkeley 
        !           495: *> \author Univ. of Colorado Denver 
        !           496: *> \author NAG Ltd. 
        !           497: *
        !           498: *> \date November 2011
        !           499: *
        !           500: *> \ingroup doubleSYcomputational
        !           501: *
        !           502: *  =====================================================================
1.1       bertrand  503:       SUBROUTINE DSYSVXX( FACT, UPLO, N, NRHS, A, LDA, AF, LDAF, IPIV,
                    504:      $                    EQUED, S, B, LDB, X, LDX, RCOND, RPVGRW, BERR,
                    505:      $                    N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
                    506:      $                    NPARAMS, PARAMS, WORK, IWORK, INFO )
                    507: *
1.5     ! bertrand  508: *  -- LAPACK computational routine (version 3.4.0) --
        !           509: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
        !           510: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
        !           511: *     November 2011
1.1       bertrand  512: *
                    513: *     .. Scalar Arguments ..
                    514:       CHARACTER          EQUED, FACT, UPLO
                    515:       INTEGER            INFO, LDA, LDAF, LDB, LDX, N, NRHS, NPARAMS,
                    516:      $                   N_ERR_BNDS
                    517:       DOUBLE PRECISION   RCOND, RPVGRW
                    518: *     ..
                    519: *     .. Array Arguments ..
                    520:       INTEGER            IPIV( * ), IWORK( * )
                    521:       DOUBLE PRECISION   A( LDA, * ), AF( LDAF, * ), B( LDB, * ),
                    522:      $                   X( LDX, * ), WORK( * )
                    523:       DOUBLE PRECISION   S( * ), PARAMS( * ), BERR( * ),
                    524:      $                   ERR_BNDS_NORM( NRHS, * ),
                    525:      $                   ERR_BNDS_COMP( NRHS, * )
                    526: *     ..
                    527: *
1.5     ! bertrand  528: *  ==================================================================
1.1       bertrand  529: *
                    530: *     .. Parameters ..
                    531:       DOUBLE PRECISION   ZERO, ONE
                    532:       PARAMETER          ( ZERO = 0.0D+0, ONE = 1.0D+0 )
                    533:       INTEGER            FINAL_NRM_ERR_I, FINAL_CMP_ERR_I, BERR_I
                    534:       INTEGER            RCOND_I, NRM_RCOND_I, NRM_ERR_I, CMP_RCOND_I
                    535:       INTEGER            CMP_ERR_I, PIV_GROWTH_I
                    536:       PARAMETER          ( FINAL_NRM_ERR_I = 1, FINAL_CMP_ERR_I = 2,
                    537:      $                   BERR_I = 3 )
                    538:       PARAMETER          ( RCOND_I = 4, NRM_RCOND_I = 5, NRM_ERR_I = 6 )
                    539:       PARAMETER          ( CMP_RCOND_I = 7, CMP_ERR_I = 8,
                    540:      $                   PIV_GROWTH_I = 9 )
                    541: *     ..
                    542: *     .. Local Scalars ..
                    543:       LOGICAL            EQUIL, NOFACT, RCEQU
                    544:       INTEGER            INFEQU, J
                    545:       DOUBLE PRECISION   AMAX, BIGNUM, SMIN, SMAX, SCOND, SMLNUM
                    546: *     ..
                    547: *     .. External Functions ..
                    548:       EXTERNAL           LSAME, DLAMCH, DLA_SYRPVGRW
                    549:       LOGICAL            LSAME
                    550:       DOUBLE PRECISION   DLAMCH, DLA_SYRPVGRW
                    551: *     ..
                    552: *     .. External Subroutines ..
                    553:       EXTERNAL           DSYCON, DSYEQUB, DSYTRF, DSYTRS,
                    554:      $                   DLACPY, DLAQSY, XERBLA, DLASCL2, DSYRFSX
                    555: *     ..
                    556: *     .. Intrinsic Functions ..
                    557:       INTRINSIC          MAX, MIN
                    558: *     ..
                    559: *     .. Executable Statements ..
                    560: *
                    561:       INFO = 0
                    562:       NOFACT = LSAME( FACT, 'N' )
                    563:       EQUIL = LSAME( FACT, 'E' )
                    564:       SMLNUM = DLAMCH( 'Safe minimum' )
                    565:       BIGNUM = ONE / SMLNUM
                    566:       IF( NOFACT .OR. EQUIL ) THEN
                    567:          EQUED = 'N'
                    568:          RCEQU = .FALSE.
                    569:       ELSE
                    570:          RCEQU = LSAME( EQUED, 'Y' )
                    571:       ENDIF
                    572: *
                    573: *     Default is failure.  If an input parameter is wrong or
                    574: *     factorization fails, make everything look horrible.  Only the
                    575: *     pivot growth is set here, the rest is initialized in DSYRFSX.
                    576: *
                    577:       RPVGRW = ZERO
                    578: *
                    579: *     Test the input parameters.  PARAMS is not tested until DSYRFSX.
                    580: *
                    581:       IF( .NOT.NOFACT .AND. .NOT.EQUIL .AND. .NOT.
                    582:      $     LSAME( FACT, 'F' ) ) THEN
                    583:          INFO = -1
                    584:       ELSE IF( .NOT.LSAME(UPLO, 'U') .AND.
                    585:      $         .NOT.LSAME(UPLO, 'L') ) THEN
                    586:          INFO = -2
                    587:       ELSE IF( N.LT.0 ) THEN
                    588:          INFO = -3
                    589:       ELSE IF( NRHS.LT.0 ) THEN
                    590:          INFO = -4
                    591:       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
                    592:          INFO = -6
                    593:       ELSE IF( LDAF.LT.MAX( 1, N ) ) THEN
                    594:          INFO = -8
                    595:       ELSE IF( LSAME( FACT, 'F' ) .AND. .NOT.
                    596:      $        ( RCEQU .OR. LSAME( EQUED, 'N' ) ) ) THEN
                    597:          INFO = -9
                    598:       ELSE
                    599:          IF ( RCEQU ) THEN
                    600:             SMIN = BIGNUM
                    601:             SMAX = ZERO
                    602:             DO 10 J = 1, N
                    603:                SMIN = MIN( SMIN, S( J ) )
                    604:                SMAX = MAX( SMAX, S( J ) )
                    605:  10         CONTINUE
                    606:             IF( SMIN.LE.ZERO ) THEN
                    607:                INFO = -10
                    608:             ELSE IF( N.GT.0 ) THEN
                    609:                SCOND = MAX( SMIN, SMLNUM ) / MIN( SMAX, BIGNUM )
                    610:             ELSE
                    611:                SCOND = ONE
                    612:             END IF
                    613:          END IF
                    614:          IF( INFO.EQ.0 ) THEN
                    615:             IF( LDB.LT.MAX( 1, N ) ) THEN
                    616:                INFO = -12
                    617:             ELSE IF( LDX.LT.MAX( 1, N ) ) THEN
                    618:                INFO = -14
                    619:             END IF
                    620:          END IF
                    621:       END IF
                    622: *
                    623:       IF( INFO.NE.0 ) THEN
                    624:          CALL XERBLA( 'DSYSVXX', -INFO )
                    625:          RETURN
                    626:       END IF
                    627: *
                    628:       IF( EQUIL ) THEN
                    629: *
                    630: *     Compute row and column scalings to equilibrate the matrix A.
                    631: *
                    632:          CALL DSYEQUB( UPLO, N, A, LDA, S, SCOND, AMAX, WORK, INFEQU )
                    633:          IF( INFEQU.EQ.0 ) THEN
                    634: *
                    635: *     Equilibrate the matrix.
                    636: *
                    637:             CALL DLAQSY( UPLO, N, A, LDA, S, SCOND, AMAX, EQUED )
                    638:             RCEQU = LSAME( EQUED, 'Y' )
                    639:          END IF
                    640:       END IF
                    641: *
                    642: *     Scale the right-hand side.
                    643: *
                    644:       IF( RCEQU ) CALL DLASCL2( N, NRHS, S, B, LDB )
                    645: *
                    646:       IF( NOFACT .OR. EQUIL ) THEN
                    647: *
                    648: *        Compute the LDL^T or UDU^T factorization of A.
                    649: *
                    650:          CALL DLACPY( UPLO, N, N, A, LDA, AF, LDAF )
                    651:          CALL DSYTRF( UPLO, N, AF, LDAF, IPIV, WORK, 5*MAX(1,N), INFO )
                    652: *
                    653: *        Return if INFO is non-zero.
                    654: *
                    655:          IF( INFO.GT.0 ) THEN
                    656: *
                    657: *           Pivot in column INFO is exactly 0
                    658: *           Compute the reciprocal pivot growth factor of the
                    659: *           leading rank-deficient INFO columns of A.
                    660: *
                    661:             IF ( N.GT.0 )
                    662:      $           RPVGRW = DLA_SYRPVGRW(UPLO, N, INFO, A, LDA, AF,
                    663:      $           LDAF, IPIV, WORK )
                    664:             RETURN
                    665:          END IF
                    666:       END IF
                    667: *
                    668: *     Compute the reciprocal pivot growth factor RPVGRW.
                    669: *
                    670:       IF ( N.GT.0 )
                    671:      $     RPVGRW = DLA_SYRPVGRW( UPLO, N, INFO, A, LDA, AF, LDAF,
                    672:      $     IPIV, WORK )
                    673: *
                    674: *     Compute the solution matrix X.
                    675: *
                    676:       CALL DLACPY( 'Full', N, NRHS, B, LDB, X, LDX )
                    677:       CALL DSYTRS( UPLO, N, NRHS, AF, LDAF, IPIV, X, LDX, INFO )
                    678: *
                    679: *     Use iterative refinement to improve the computed solution and
                    680: *     compute error bounds and backward error estimates for it.
                    681: *
                    682:       CALL DSYRFSX( UPLO, EQUED, N, NRHS, A, LDA, AF, LDAF, IPIV,
                    683:      $     S, B, LDB, X, LDX, RCOND, BERR, N_ERR_BNDS, ERR_BNDS_NORM,
                    684:      $     ERR_BNDS_COMP, NPARAMS, PARAMS, WORK, IWORK, INFO )
                    685: *
                    686: *     Scale solutions.
                    687: *
                    688:       IF ( RCEQU ) THEN
                    689:          CALL DLASCL2 ( N, NRHS, S, X, LDX )
                    690:       END IF
                    691: *
                    692:       RETURN
                    693: *
                    694: *     End of DSYSVXX
                    695: *
                    696:       END

CVSweb interface <joel.bertrand@systella.fr>