File:  [local] / rpl / lapack / lapack / dgetrs.f
Revision 1.12: download - view: text, annotated - select for diffs - revision graph
Fri Dec 14 14:22:30 2012 UTC (11 years, 5 months ago) by bertrand
Branches: MAIN
CVS tags: rpl-4_1_16, rpl-4_1_15, rpl-4_1_14, rpl-4_1_13, rpl-4_1_12, rpl-4_1_11, HEAD
Mise à jour de lapack.

    1: *> \brief \b DGETRS
    2: *
    3: *  =========== DOCUMENTATION ===========
    4: *
    5: * Online html documentation available at 
    6: *            http://www.netlib.org/lapack/explore-html/ 
    7: *
    8: *> \htmlonly
    9: *> Download DGETRS + dependencies 
   10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgetrs.f"> 
   11: *> [TGZ]</a> 
   12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgetrs.f"> 
   13: *> [ZIP]</a> 
   14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgetrs.f"> 
   15: *> [TXT]</a>
   16: *> \endhtmlonly 
   17: *
   18: *  Definition:
   19: *  ===========
   20: *
   21: *       SUBROUTINE DGETRS( TRANS, N, NRHS, A, LDA, IPIV, B, LDB, INFO )
   22:    23: *       .. Scalar Arguments ..
   24: *       CHARACTER          TRANS
   25: *       INTEGER            INFO, LDA, LDB, N, NRHS
   26: *       ..
   27: *       .. Array Arguments ..
   28: *       INTEGER            IPIV( * )
   29: *       DOUBLE PRECISION   A( LDA, * ), B( LDB, * )
   30: *       ..
   31: *  
   32: *
   33: *> \par Purpose:
   34: *  =============
   35: *>
   36: *> \verbatim
   37: *>
   38: *> DGETRS solves a system of linear equations
   39: *>    A * X = B  or  A**T * X = B
   40: *> with a general N-by-N matrix A using the LU factorization computed
   41: *> by DGETRF.
   42: *> \endverbatim
   43: *
   44: *  Arguments:
   45: *  ==========
   46: *
   47: *> \param[in] TRANS
   48: *> \verbatim
   49: *>          TRANS is CHARACTER*1
   50: *>          Specifies the form of the system of equations:
   51: *>          = 'N':  A * X = B  (No transpose)
   52: *>          = 'T':  A**T* X = B  (Transpose)
   53: *>          = 'C':  A**T* X = B  (Conjugate transpose = Transpose)
   54: *> \endverbatim
   55: *>
   56: *> \param[in] N
   57: *> \verbatim
   58: *>          N is INTEGER
   59: *>          The order of the matrix A.  N >= 0.
   60: *> \endverbatim
   61: *>
   62: *> \param[in] NRHS
   63: *> \verbatim
   64: *>          NRHS is INTEGER
   65: *>          The number of right hand sides, i.e., the number of columns
   66: *>          of the matrix B.  NRHS >= 0.
   67: *> \endverbatim
   68: *>
   69: *> \param[in] A
   70: *> \verbatim
   71: *>          A is DOUBLE PRECISION array, dimension (LDA,N)
   72: *>          The factors L and U from the factorization A = P*L*U
   73: *>          as computed by DGETRF.
   74: *> \endverbatim
   75: *>
   76: *> \param[in] LDA
   77: *> \verbatim
   78: *>          LDA is INTEGER
   79: *>          The leading dimension of the array A.  LDA >= max(1,N).
   80: *> \endverbatim
   81: *>
   82: *> \param[in] IPIV
   83: *> \verbatim
   84: *>          IPIV is INTEGER array, dimension (N)
   85: *>          The pivot indices from DGETRF; for 1<=i<=N, row i of the
   86: *>          matrix was interchanged with row IPIV(i).
   87: *> \endverbatim
   88: *>
   89: *> \param[in,out] B
   90: *> \verbatim
   91: *>          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
   92: *>          On entry, the right hand side matrix B.
   93: *>          On exit, the solution matrix X.
   94: *> \endverbatim
   95: *>
   96: *> \param[in] LDB
   97: *> \verbatim
   98: *>          LDB is INTEGER
   99: *>          The leading dimension of the array B.  LDB >= max(1,N).
  100: *> \endverbatim
  101: *>
  102: *> \param[out] INFO
  103: *> \verbatim
  104: *>          INFO is INTEGER
  105: *>          = 0:  successful exit
  106: *>          < 0:  if INFO = -i, the i-th argument had an illegal value
  107: *> \endverbatim
  108: *
  109: *  Authors:
  110: *  ========
  111: *
  112: *> \author Univ. of Tennessee 
  113: *> \author Univ. of California Berkeley 
  114: *> \author Univ. of Colorado Denver 
  115: *> \author NAG Ltd. 
  116: *
  117: *> \date November 2011
  118: *
  119: *> \ingroup doubleGEcomputational
  120: *
  121: *  =====================================================================
  122:       SUBROUTINE DGETRS( TRANS, N, NRHS, A, LDA, IPIV, B, LDB, INFO )
  123: *
  124: *  -- LAPACK computational routine (version 3.4.0) --
  125: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
  126: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  127: *     November 2011
  128: *
  129: *     .. Scalar Arguments ..
  130:       CHARACTER          TRANS
  131:       INTEGER            INFO, LDA, LDB, N, NRHS
  132: *     ..
  133: *     .. Array Arguments ..
  134:       INTEGER            IPIV( * )
  135:       DOUBLE PRECISION   A( LDA, * ), B( LDB, * )
  136: *     ..
  137: *
  138: *  =====================================================================
  139: *
  140: *     .. Parameters ..
  141:       DOUBLE PRECISION   ONE
  142:       PARAMETER          ( ONE = 1.0D+0 )
  143: *     ..
  144: *     .. Local Scalars ..
  145:       LOGICAL            NOTRAN
  146: *     ..
  147: *     .. External Functions ..
  148:       LOGICAL            LSAME
  149:       EXTERNAL           LSAME
  150: *     ..
  151: *     .. External Subroutines ..
  152:       EXTERNAL           DLASWP, DTRSM, XERBLA
  153: *     ..
  154: *     .. Intrinsic Functions ..
  155:       INTRINSIC          MAX
  156: *     ..
  157: *     .. Executable Statements ..
  158: *
  159: *     Test the input parameters.
  160: *
  161:       INFO = 0
  162:       NOTRAN = LSAME( TRANS, 'N' )
  163:       IF( .NOT.NOTRAN .AND. .NOT.LSAME( TRANS, 'T' ) .AND. .NOT.
  164:      $    LSAME( TRANS, 'C' ) ) THEN
  165:          INFO = -1
  166:       ELSE IF( N.LT.0 ) THEN
  167:          INFO = -2
  168:       ELSE IF( NRHS.LT.0 ) THEN
  169:          INFO = -3
  170:       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
  171:          INFO = -5
  172:       ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
  173:          INFO = -8
  174:       END IF
  175:       IF( INFO.NE.0 ) THEN
  176:          CALL XERBLA( 'DGETRS', -INFO )
  177:          RETURN
  178:       END IF
  179: *
  180: *     Quick return if possible
  181: *
  182:       IF( N.EQ.0 .OR. NRHS.EQ.0 )
  183:      $   RETURN
  184: *
  185:       IF( NOTRAN ) THEN
  186: *
  187: *        Solve A * X = B.
  188: *
  189: *        Apply row interchanges to the right hand sides.
  190: *
  191:          CALL DLASWP( NRHS, B, LDB, 1, N, IPIV, 1 )
  192: *
  193: *        Solve L*X = B, overwriting B with X.
  194: *
  195:          CALL DTRSM( 'Left', 'Lower', 'No transpose', 'Unit', N, NRHS,
  196:      $               ONE, A, LDA, B, LDB )
  197: *
  198: *        Solve U*X = B, overwriting B with X.
  199: *
  200:          CALL DTRSM( 'Left', 'Upper', 'No transpose', 'Non-unit', N,
  201:      $               NRHS, ONE, A, LDA, B, LDB )
  202:       ELSE
  203: *
  204: *        Solve A**T * X = B.
  205: *
  206: *        Solve U**T *X = B, overwriting B with X.
  207: *
  208:          CALL DTRSM( 'Left', 'Upper', 'Transpose', 'Non-unit', N, NRHS,
  209:      $               ONE, A, LDA, B, LDB )
  210: *
  211: *        Solve L**T *X = B, overwriting B with X.
  212: *
  213:          CALL DTRSM( 'Left', 'Lower', 'Transpose', 'Unit', N, NRHS, ONE,
  214:      $               A, LDA, B, LDB )
  215: *
  216: *        Apply row interchanges to the solution vectors.
  217: *
  218:          CALL DLASWP( NRHS, B, LDB, 1, N, IPIV, -1 )
  219:       END IF
  220: *
  221:       RETURN
  222: *
  223: *     End of DGETRS
  224: *
  225:       END

CVSweb interface <joel.bertrand@systella.fr>