File:  [local] / rpl / lapack / lapack / dgesv.f
Revision 1.17: download - view: text, annotated - select for diffs - revision graph
Mon Aug 7 08:38:50 2023 UTC (8 months, 3 weeks ago) by bertrand
Branches: MAIN
CVS tags: rpl-4_1_35, rpl-4_1_34, HEAD
Première mise à jour de lapack et blas.

    1: *> \brief <b> DGESV computes the solution to system of linear equations A * X = B for GE matrices</b>
    2: *
    3: *  =========== DOCUMENTATION ===========
    4: *
    5: * Online html documentation available at
    6: *            http://www.netlib.org/lapack/explore-html/
    7: *
    8: *> \htmlonly
    9: *> Download DGESV + dependencies
   10: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgesv.f">
   11: *> [TGZ]</a>
   12: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgesv.f">
   13: *> [ZIP]</a>
   14: *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgesv.f">
   15: *> [TXT]</a>
   16: *> \endhtmlonly
   17: *
   18: *  Definition:
   19: *  ===========
   20: *
   21: *       SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO )
   22: *
   23: *       .. Scalar Arguments ..
   24: *       INTEGER            INFO, LDA, LDB, N, NRHS
   25: *       ..
   26: *       .. Array Arguments ..
   27: *       INTEGER            IPIV( * )
   28: *       DOUBLE PRECISION   A( LDA, * ), B( LDB, * )
   29: *       ..
   30: *
   31: *
   32: *> \par Purpose:
   33: *  =============
   34: *>
   35: *> \verbatim
   36: *>
   37: *> DGESV computes the solution to a real system of linear equations
   38: *>    A * X = B,
   39: *> where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
   40: *>
   41: *> The LU decomposition with partial pivoting and row interchanges is
   42: *> used to factor A as
   43: *>    A = P * L * U,
   44: *> where P is a permutation matrix, L is unit lower triangular, and U is
   45: *> upper triangular.  The factored form of A is then used to solve the
   46: *> system of equations A * X = B.
   47: *> \endverbatim
   48: *
   49: *  Arguments:
   50: *  ==========
   51: *
   52: *> \param[in] N
   53: *> \verbatim
   54: *>          N is INTEGER
   55: *>          The number of linear equations, i.e., the order of the
   56: *>          matrix A.  N >= 0.
   57: *> \endverbatim
   58: *>
   59: *> \param[in] NRHS
   60: *> \verbatim
   61: *>          NRHS is INTEGER
   62: *>          The number of right hand sides, i.e., the number of columns
   63: *>          of the matrix B.  NRHS >= 0.
   64: *> \endverbatim
   65: *>
   66: *> \param[in,out] A
   67: *> \verbatim
   68: *>          A is DOUBLE PRECISION array, dimension (LDA,N)
   69: *>          On entry, the N-by-N coefficient matrix A.
   70: *>          On exit, the factors L and U from the factorization
   71: *>          A = P*L*U; the unit diagonal elements of L are not stored.
   72: *> \endverbatim
   73: *>
   74: *> \param[in] LDA
   75: *> \verbatim
   76: *>          LDA is INTEGER
   77: *>          The leading dimension of the array A.  LDA >= max(1,N).
   78: *> \endverbatim
   79: *>
   80: *> \param[out] IPIV
   81: *> \verbatim
   82: *>          IPIV is INTEGER array, dimension (N)
   83: *>          The pivot indices that define the permutation matrix P;
   84: *>          row i of the matrix was interchanged with row IPIV(i).
   85: *> \endverbatim
   86: *>
   87: *> \param[in,out] B
   88: *> \verbatim
   89: *>          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
   90: *>          On entry, the N-by-NRHS matrix of right hand side matrix B.
   91: *>          On exit, if INFO = 0, the N-by-NRHS solution matrix X.
   92: *> \endverbatim
   93: *>
   94: *> \param[in] LDB
   95: *> \verbatim
   96: *>          LDB is INTEGER
   97: *>          The leading dimension of the array B.  LDB >= max(1,N).
   98: *> \endverbatim
   99: *>
  100: *> \param[out] INFO
  101: *> \verbatim
  102: *>          INFO is INTEGER
  103: *>          = 0:  successful exit
  104: *>          < 0:  if INFO = -i, the i-th argument had an illegal value
  105: *>          > 0:  if INFO = i, U(i,i) is exactly zero.  The factorization
  106: *>                has been completed, but the factor U is exactly
  107: *>                singular, so the solution could not be computed.
  108: *> \endverbatim
  109: *
  110: *  Authors:
  111: *  ========
  112: *
  113: *> \author Univ. of Tennessee
  114: *> \author Univ. of California Berkeley
  115: *> \author Univ. of Colorado Denver
  116: *> \author NAG Ltd.
  117: *
  118: *> \ingroup doubleGEsolve
  119: *
  120: *  =====================================================================
  121:       SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO )
  122: *
  123: *  -- LAPACK driver routine --
  124: *  -- LAPACK is a software package provided by Univ. of Tennessee,    --
  125: *  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
  126: *
  127: *     .. Scalar Arguments ..
  128:       INTEGER            INFO, LDA, LDB, N, NRHS
  129: *     ..
  130: *     .. Array Arguments ..
  131:       INTEGER            IPIV( * )
  132:       DOUBLE PRECISION   A( LDA, * ), B( LDB, * )
  133: *     ..
  134: *
  135: *  =====================================================================
  136: *
  137: *     .. External Subroutines ..
  138:       EXTERNAL           DGETRF, DGETRS, XERBLA
  139: *     ..
  140: *     .. Intrinsic Functions ..
  141:       INTRINSIC          MAX
  142: *     ..
  143: *     .. Executable Statements ..
  144: *
  145: *     Test the input parameters.
  146: *
  147:       INFO = 0
  148:       IF( N.LT.0 ) THEN
  149:          INFO = -1
  150:       ELSE IF( NRHS.LT.0 ) THEN
  151:          INFO = -2
  152:       ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
  153:          INFO = -4
  154:       ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
  155:          INFO = -7
  156:       END IF
  157:       IF( INFO.NE.0 ) THEN
  158:          CALL XERBLA( 'DGESV ', -INFO )
  159:          RETURN
  160:       END IF
  161: *
  162: *     Compute the LU factorization of A.
  163: *
  164:       CALL DGETRF( N, N, A, LDA, IPIV, INFO )
  165:       IF( INFO.EQ.0 ) THEN
  166: *
  167: *        Solve the system A*X = B, overwriting B with X.
  168: *
  169:          CALL DGETRS( 'No transpose', N, NRHS, A, LDA, IPIV, B, LDB,
  170:      $                INFO )
  171:       END IF
  172:       RETURN
  173: *
  174: *     End of DGESV
  175: *
  176:       END

CVSweb interface <joel.bertrand@systella.fr>