2017-10-25 01:22:20 +08:00
|
|
|
*> \brief <b> DGESV computes the solution to system of linear equations A * X = B for GE matrices</b>
|
|
|
|
*
|
|
|
|
* =========== DOCUMENTATION ===========
|
|
|
|
*
|
2018-05-19 05:17:13 +08:00
|
|
|
* Online html documentation available at
|
|
|
|
* http://www.netlib.org/lapack/explore-html/
|
2017-10-25 01:22:20 +08:00
|
|
|
*
|
|
|
|
*> \htmlonly
|
2018-05-19 05:17:13 +08:00
|
|
|
*> Download DGESV + dependencies
|
|
|
|
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgesv.f">
|
|
|
|
*> [TGZ]</a>
|
|
|
|
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgesv.f">
|
|
|
|
*> [ZIP]</a>
|
|
|
|
*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgesv.f">
|
2017-10-25 01:22:20 +08:00
|
|
|
*> [TXT]</a>
|
|
|
|
*> \endhtmlonly
|
|
|
|
*
|
|
|
|
* Definition:
|
|
|
|
* ===========
|
|
|
|
*
|
|
|
|
* SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO )
|
2018-05-19 05:17:13 +08:00
|
|
|
*
|
2017-10-25 01:22:20 +08:00
|
|
|
* .. Scalar Arguments ..
|
|
|
|
* INTEGER INFO, LDA, LDB, N, NRHS
|
|
|
|
* ..
|
|
|
|
* .. Array Arguments ..
|
|
|
|
* INTEGER IPIV( * )
|
|
|
|
* DOUBLE PRECISION A( LDA, * ), B( LDB, * )
|
|
|
|
* ..
|
2018-05-19 05:17:13 +08:00
|
|
|
*
|
2017-10-25 01:22:20 +08:00
|
|
|
*
|
|
|
|
*> \par Purpose:
|
|
|
|
* =============
|
|
|
|
*>
|
|
|
|
*> \verbatim
|
|
|
|
*>
|
|
|
|
*> DGESV computes the solution to a real system of linear equations
|
|
|
|
*> A * X = B,
|
|
|
|
*> where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
|
|
|
|
*>
|
|
|
|
*> The LU decomposition with partial pivoting and row interchanges is
|
|
|
|
*> used to factor A as
|
|
|
|
*> A = P * L * U,
|
|
|
|
*> where P is a permutation matrix, L is unit lower triangular, and U is
|
|
|
|
*> upper triangular. The factored form of A is then used to solve the
|
|
|
|
*> system of equations A * X = B.
|
|
|
|
*> \endverbatim
|
|
|
|
*
|
|
|
|
* Arguments:
|
|
|
|
* ==========
|
|
|
|
*
|
|
|
|
*> \param[in] N
|
|
|
|
*> \verbatim
|
|
|
|
*> N is INTEGER
|
|
|
|
*> The number of linear equations, i.e., the order of the
|
|
|
|
*> matrix A. N >= 0.
|
|
|
|
*> \endverbatim
|
|
|
|
*>
|
|
|
|
*> \param[in] NRHS
|
|
|
|
*> \verbatim
|
|
|
|
*> NRHS is INTEGER
|
|
|
|
*> The number of right hand sides, i.e., the number of columns
|
|
|
|
*> of the matrix B. NRHS >= 0.
|
|
|
|
*> \endverbatim
|
|
|
|
*>
|
|
|
|
*> \param[in,out] A
|
|
|
|
*> \verbatim
|
|
|
|
*> A is DOUBLE PRECISION array, dimension (LDA,N)
|
|
|
|
*> On entry, the N-by-N coefficient matrix A.
|
|
|
|
*> On exit, the factors L and U from the factorization
|
|
|
|
*> A = P*L*U; the unit diagonal elements of L are not stored.
|
|
|
|
*> \endverbatim
|
|
|
|
*>
|
|
|
|
*> \param[in] LDA
|
|
|
|
*> \verbatim
|
|
|
|
*> LDA is INTEGER
|
|
|
|
*> The leading dimension of the array A. LDA >= max(1,N).
|
|
|
|
*> \endverbatim
|
|
|
|
*>
|
|
|
|
*> \param[out] IPIV
|
|
|
|
*> \verbatim
|
|
|
|
*> IPIV is INTEGER array, dimension (N)
|
|
|
|
*> The pivot indices that define the permutation matrix P;
|
|
|
|
*> row i of the matrix was interchanged with row IPIV(i).
|
|
|
|
*> \endverbatim
|
|
|
|
*>
|
|
|
|
*> \param[in,out] B
|
|
|
|
*> \verbatim
|
|
|
|
*> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
|
|
|
|
*> On entry, the N-by-NRHS matrix of right hand side matrix B.
|
|
|
|
*> On exit, if INFO = 0, the N-by-NRHS solution matrix X.
|
|
|
|
*> \endverbatim
|
|
|
|
*>
|
|
|
|
*> \param[in] LDB
|
|
|
|
*> \verbatim
|
|
|
|
*> LDB is INTEGER
|
|
|
|
*> The leading dimension of the array B. LDB >= max(1,N).
|
|
|
|
*> \endverbatim
|
|
|
|
*>
|
|
|
|
*> \param[out] INFO
|
|
|
|
*> \verbatim
|
|
|
|
*> INFO is INTEGER
|
|
|
|
*> = 0: successful exit
|
|
|
|
*> < 0: if INFO = -i, the i-th argument had an illegal value
|
|
|
|
*> > 0: if INFO = i, U(i,i) is exactly zero. The factorization
|
|
|
|
*> has been completed, but the factor U is exactly
|
|
|
|
*> singular, so the solution could not be computed.
|
|
|
|
*> \endverbatim
|
|
|
|
*
|
|
|
|
* Authors:
|
|
|
|
* ========
|
|
|
|
*
|
2018-05-19 05:17:13 +08:00
|
|
|
*> \author Univ. of Tennessee
|
|
|
|
*> \author Univ. of California Berkeley
|
|
|
|
*> \author Univ. of Colorado Denver
|
|
|
|
*> \author NAG Ltd.
|
2017-10-25 01:22:20 +08:00
|
|
|
*
|
2018-05-19 05:17:13 +08:00
|
|
|
*> \date December 2016
|
2017-10-25 01:22:20 +08:00
|
|
|
*
|
|
|
|
*> \ingroup doubleGEsolve
|
|
|
|
*
|
|
|
|
* =====================================================================
|
|
|
|
SUBROUTINE DGESV( N, NRHS, A, LDA, IPIV, B, LDB, INFO )
|
|
|
|
*
|
2018-05-19 05:17:13 +08:00
|
|
|
* -- LAPACK driver routine (version 3.7.0) --
|
2017-10-25 01:22:20 +08:00
|
|
|
* -- LAPACK is a software package provided by Univ. of Tennessee, --
|
|
|
|
* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
|
2018-05-19 05:17:13 +08:00
|
|
|
* December 2016
|
2017-10-25 01:22:20 +08:00
|
|
|
*
|
|
|
|
* .. Scalar Arguments ..
|
|
|
|
INTEGER INFO, LDA, LDB, N, NRHS
|
|
|
|
* ..
|
|
|
|
* .. Array Arguments ..
|
|
|
|
INTEGER IPIV( * )
|
|
|
|
DOUBLE PRECISION A( LDA, * ), B( LDB, * )
|
|
|
|
* ..
|
|
|
|
*
|
|
|
|
* =====================================================================
|
|
|
|
*
|
|
|
|
* .. External Subroutines ..
|
|
|
|
EXTERNAL DGETRF, DGETRS, XERBLA
|
|
|
|
* ..
|
|
|
|
* .. Intrinsic Functions ..
|
|
|
|
INTRINSIC MAX
|
|
|
|
* ..
|
|
|
|
* .. Executable Statements ..
|
|
|
|
*
|
|
|
|
* Test the input parameters.
|
|
|
|
*
|
|
|
|
INFO = 0
|
|
|
|
IF( N.LT.0 ) THEN
|
|
|
|
INFO = -1
|
|
|
|
ELSE IF( NRHS.LT.0 ) THEN
|
|
|
|
INFO = -2
|
|
|
|
ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
|
|
|
|
INFO = -4
|
|
|
|
ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
|
|
|
|
INFO = -7
|
|
|
|
END IF
|
|
|
|
IF( INFO.NE.0 ) THEN
|
|
|
|
CALL XERBLA( 'DGESV ', -INFO )
|
|
|
|
RETURN
|
|
|
|
END IF
|
|
|
|
*
|
|
|
|
* Compute the LU factorization of A.
|
|
|
|
*
|
|
|
|
CALL DGETRF( N, N, A, LDA, IPIV, INFO )
|
|
|
|
IF( INFO.EQ.0 ) THEN
|
|
|
|
*
|
|
|
|
* Solve the system A*X = B, overwriting B with X.
|
|
|
|
*
|
|
|
|
CALL DGETRS( 'No transpose', N, NRHS, A, LDA, IPIV, B, LDB,
|
|
|
|
$ INFO )
|
|
|
|
END IF
|
|
|
|
RETURN
|
|
|
|
*
|
|
|
|
* End of DGESV
|
|
|
|
*
|
|
|
|
END
|