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: *> \ingroup doubleGEcomputational
118: *
119: * =====================================================================
120: SUBROUTINE DGETRS( TRANS, N, NRHS, A, LDA, IPIV, B, LDB, INFO )
121: *
122: * -- LAPACK computational routine --
123: * -- LAPACK is a software package provided by Univ. of Tennessee, --
124: * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
125: *
126: * .. Scalar Arguments ..
127: CHARACTER TRANS
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: * .. Parameters ..
138: DOUBLE PRECISION ONE
139: PARAMETER ( ONE = 1.0D+0 )
140: * ..
141: * .. Local Scalars ..
142: LOGICAL NOTRAN
143: * ..
144: * .. External Functions ..
145: LOGICAL LSAME
146: EXTERNAL LSAME
147: * ..
148: * .. External Subroutines ..
149: EXTERNAL DLASWP, DTRSM, XERBLA
150: * ..
151: * .. Intrinsic Functions ..
152: INTRINSIC MAX
153: * ..
154: * .. Executable Statements ..
155: *
156: * Test the input parameters.
157: *
158: INFO = 0
159: NOTRAN = LSAME( TRANS, 'N' )
160: IF( .NOT.NOTRAN .AND. .NOT.LSAME( TRANS, 'T' ) .AND. .NOT.
161: $ LSAME( TRANS, 'C' ) ) THEN
162: INFO = -1
163: ELSE IF( N.LT.0 ) THEN
164: INFO = -2
165: ELSE IF( NRHS.LT.0 ) THEN
166: INFO = -3
167: ELSE IF( LDA.LT.MAX( 1, N ) ) THEN
168: INFO = -5
169: ELSE IF( LDB.LT.MAX( 1, N ) ) THEN
170: INFO = -8
171: END IF
172: IF( INFO.NE.0 ) THEN
173: CALL XERBLA( 'DGETRS', -INFO )
174: RETURN
175: END IF
176: *
177: * Quick return if possible
178: *
179: IF( N.EQ.0 .OR. NRHS.EQ.0 )
180: $ RETURN
181: *
182: IF( NOTRAN ) THEN
183: *
184: * Solve A * X = B.
185: *
186: * Apply row interchanges to the right hand sides.
187: *
188: CALL DLASWP( NRHS, B, LDB, 1, N, IPIV, 1 )
189: *
190: * Solve L*X = B, overwriting B with X.
191: *
192: CALL DTRSM( 'Left', 'Lower', 'No transpose', 'Unit', N, NRHS,
193: $ ONE, A, LDA, B, LDB )
194: *
195: * Solve U*X = B, overwriting B with X.
196: *
197: CALL DTRSM( 'Left', 'Upper', 'No transpose', 'Non-unit', N,
198: $ NRHS, ONE, A, LDA, B, LDB )
199: ELSE
200: *
201: * Solve A**T * X = B.
202: *
203: * Solve U**T *X = B, overwriting B with X.
204: *
205: CALL DTRSM( 'Left', 'Upper', 'Transpose', 'Non-unit', N, NRHS,
206: $ ONE, A, LDA, B, LDB )
207: *
208: * Solve L**T *X = B, overwriting B with X.
209: *
210: CALL DTRSM( 'Left', 'Lower', 'Transpose', 'Unit', N, NRHS, ONE,
211: $ A, LDA, B, LDB )
212: *
213: * Apply row interchanges to the solution vectors.
214: *
215: CALL DLASWP( NRHS, B, LDB, 1, N, IPIV, -1 )
216: END IF
217: *
218: RETURN
219: *
220: * End of DGETRS
221: *
222: END
CVSweb interface <joel.bertrand@systella.fr>