LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dgtsv.f
Go to the documentation of this file.
1*> \brief <b> DGTSV computes the solution to system of linear equations A * X = B for GT matrices </b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download DGTSV + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgtsv.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgtsv.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgtsv.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE DGTSV( N, NRHS, DL, D, DU, B, LDB, INFO )
20*
21* .. Scalar Arguments ..
22* INTEGER INFO, LDB, N, NRHS
23* ..
24* .. Array Arguments ..
25* DOUBLE PRECISION B( LDB, * ), D( * ), DL( * ), DU( * )
26* ..
27*
28*
29*> \par Purpose:
30* =============
31*>
32*> \verbatim
33*>
34*> DGTSV solves the equation
35*>
36*> A*X = B,
37*>
38*> where A is an n by n tridiagonal matrix, by Gaussian elimination with
39*> partial pivoting.
40*>
41*> Note that the equation A**T*X = B may be solved by interchanging the
42*> order of the arguments DU and DL.
43*> \endverbatim
44*
45* Arguments:
46* ==========
47*
48*> \param[in] N
49*> \verbatim
50*> N is INTEGER
51*> The order of the matrix A. N >= 0.
52*> \endverbatim
53*>
54*> \param[in] NRHS
55*> \verbatim
56*> NRHS is INTEGER
57*> The number of right hand sides, i.e., the number of columns
58*> of the matrix B. NRHS >= 0.
59*> \endverbatim
60*>
61*> \param[in,out] DL
62*> \verbatim
63*> DL is DOUBLE PRECISION array, dimension (N-1)
64*> On entry, DL must contain the (n-1) sub-diagonal elements of
65*> A.
66*>
67*> On exit, DL is overwritten by the (n-2) elements of the
68*> second super-diagonal of the upper triangular matrix U from
69*> the LU factorization of A, in DL(1), ..., DL(n-2).
70*> \endverbatim
71*>
72*> \param[in,out] D
73*> \verbatim
74*> D is DOUBLE PRECISION array, dimension (N)
75*> On entry, D must contain the diagonal elements of A.
76*>
77*> On exit, D is overwritten by the n diagonal elements of U.
78*> \endverbatim
79*>
80*> \param[in,out] DU
81*> \verbatim
82*> DU is DOUBLE PRECISION array, dimension (N-1)
83*> On entry, DU must contain the (n-1) super-diagonal elements
84*> of A.
85*>
86*> On exit, DU is overwritten by the (n-1) elements of the first
87*> super-diagonal of U.
88*> \endverbatim
89*>
90*> \param[in,out] B
91*> \verbatim
92*> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
93*> On entry, the N by NRHS matrix of right hand side matrix B.
94*> On exit, if INFO = 0, the N by NRHS solution matrix X.
95*> \endverbatim
96*>
97*> \param[in] LDB
98*> \verbatim
99*> LDB is INTEGER
100*> The leading dimension of the array B. LDB >= max(1,N).
101*> \endverbatim
102*>
103*> \param[out] INFO
104*> \verbatim
105*> INFO is INTEGER
106*> = 0: successful exit
107*> < 0: if INFO = -i, the i-th argument had an illegal value
108*> > 0: if INFO = i, U(i,i) is exactly zero, and the solution
109*> has not been computed. The factorization has not been
110*> completed unless i = N.
111*> \endverbatim
112*
113* Authors:
114* ========
115*
116*> \author Univ. of Tennessee
117*> \author Univ. of California Berkeley
118*> \author Univ. of Colorado Denver
119*> \author NAG Ltd.
120*
121*> \ingroup gtsv
122*
123* =====================================================================
124 SUBROUTINE dgtsv( N, NRHS, DL, D, DU, B, LDB, INFO )
125*
126* -- LAPACK driver routine --
127* -- LAPACK is a software package provided by Univ. of Tennessee, --
128* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
129*
130* .. Scalar Arguments ..
131 INTEGER INFO, LDB, N, NRHS
132* ..
133* .. Array Arguments ..
134 DOUBLE PRECISION B( LDB, * ), D( * ), DL( * ), DU( * )
135* ..
136*
137* =====================================================================
138*
139* .. Parameters ..
140 DOUBLE PRECISION ZERO
141 parameter( zero = 0.0d+0 )
142* ..
143* .. Local Scalars ..
144 INTEGER I, J
145 DOUBLE PRECISION FACT, TEMP
146* ..
147* .. Intrinsic Functions ..
148 INTRINSIC abs, max
149* ..
150* .. External Subroutines ..
151 EXTERNAL xerbla
152* ..
153* .. Executable Statements ..
154*
155 info = 0
156 IF( n.LT.0 ) THEN
157 info = -1
158 ELSE IF( nrhs.LT.0 ) THEN
159 info = -2
160 ELSE IF( ldb.LT.max( 1, n ) ) THEN
161 info = -7
162 END IF
163 IF( info.NE.0 ) THEN
164 CALL xerbla( 'DGTSV ', -info )
165 RETURN
166 END IF
167*
168 IF( n.EQ.0 )
169 $ RETURN
170*
171 IF( nrhs.EQ.1 ) THEN
172 DO 10 i = 1, n - 2
173 IF( abs( d( i ) ).GE.abs( dl( i ) ) ) THEN
174*
175* No row interchange required
176*
177 IF( d( i ).NE.zero ) THEN
178 fact = dl( i ) / d( i )
179 d( i+1 ) = d( i+1 ) - fact*du( i )
180 b( i+1, 1 ) = b( i+1, 1 ) - fact*b( i, 1 )
181 ELSE
182 info = i
183 RETURN
184 END IF
185 dl( i ) = zero
186 ELSE
187*
188* Interchange rows I and I+1
189*
190 fact = d( i ) / dl( i )
191 d( i ) = dl( i )
192 temp = d( i+1 )
193 d( i+1 ) = du( i ) - fact*temp
194 dl( i ) = du( i+1 )
195 du( i+1 ) = -fact*dl( i )
196 du( i ) = temp
197 temp = b( i, 1 )
198 b( i, 1 ) = b( i+1, 1 )
199 b( i+1, 1 ) = temp - fact*b( i+1, 1 )
200 END IF
201 10 CONTINUE
202 IF( n.GT.1 ) THEN
203 i = n - 1
204 IF( abs( d( i ) ).GE.abs( dl( i ) ) ) THEN
205 IF( d( i ).NE.zero ) THEN
206 fact = dl( i ) / d( i )
207 d( i+1 ) = d( i+1 ) - fact*du( i )
208 b( i+1, 1 ) = b( i+1, 1 ) - fact*b( i, 1 )
209 ELSE
210 info = i
211 RETURN
212 END IF
213 ELSE
214 fact = d( i ) / dl( i )
215 d( i ) = dl( i )
216 temp = d( i+1 )
217 d( i+1 ) = du( i ) - fact*temp
218 du( i ) = temp
219 temp = b( i, 1 )
220 b( i, 1 ) = b( i+1, 1 )
221 b( i+1, 1 ) = temp - fact*b( i+1, 1 )
222 END IF
223 END IF
224 IF( d( n ).EQ.zero ) THEN
225 info = n
226 RETURN
227 END IF
228 ELSE
229 DO 40 i = 1, n - 2
230 IF( abs( d( i ) ).GE.abs( dl( i ) ) ) THEN
231*
232* No row interchange required
233*
234 IF( d( i ).NE.zero ) THEN
235 fact = dl( i ) / d( i )
236 d( i+1 ) = d( i+1 ) - fact*du( i )
237 DO 20 j = 1, nrhs
238 b( i+1, j ) = b( i+1, j ) - fact*b( i, j )
239 20 CONTINUE
240 ELSE
241 info = i
242 RETURN
243 END IF
244 dl( i ) = zero
245 ELSE
246*
247* Interchange rows I and I+1
248*
249 fact = d( i ) / dl( i )
250 d( i ) = dl( i )
251 temp = d( i+1 )
252 d( i+1 ) = du( i ) - fact*temp
253 dl( i ) = du( i+1 )
254 du( i+1 ) = -fact*dl( i )
255 du( i ) = temp
256 DO 30 j = 1, nrhs
257 temp = b( i, j )
258 b( i, j ) = b( i+1, j )
259 b( i+1, j ) = temp - fact*b( i+1, j )
260 30 CONTINUE
261 END IF
262 40 CONTINUE
263 IF( n.GT.1 ) THEN
264 i = n - 1
265 IF( abs( d( i ) ).GE.abs( dl( i ) ) ) THEN
266 IF( d( i ).NE.zero ) THEN
267 fact = dl( i ) / d( i )
268 d( i+1 ) = d( i+1 ) - fact*du( i )
269 DO 50 j = 1, nrhs
270 b( i+1, j ) = b( i+1, j ) - fact*b( i, j )
271 50 CONTINUE
272 ELSE
273 info = i
274 RETURN
275 END IF
276 ELSE
277 fact = d( i ) / dl( i )
278 d( i ) = dl( i )
279 temp = d( i+1 )
280 d( i+1 ) = du( i ) - fact*temp
281 du( i ) = temp
282 DO 60 j = 1, nrhs
283 temp = b( i, j )
284 b( i, j ) = b( i+1, j )
285 b( i+1, j ) = temp - fact*b( i+1, j )
286 60 CONTINUE
287 END IF
288 END IF
289 IF( d( n ).EQ.zero ) THEN
290 info = n
291 RETURN
292 END IF
293 END IF
294*
295* Back solve with the matrix U from the factorization.
296*
297 IF( nrhs.LE.2 ) THEN
298 j = 1
299 70 CONTINUE
300 b( n, j ) = b( n, j ) / d( n )
301 IF( n.GT.1 )
302 $ b( n-1, j ) = ( b( n-1, j )-du( n-1 )*b( n, j ) ) / d( n-1 )
303 DO 80 i = n - 2, 1, -1
304 b( i, j ) = ( b( i, j )-du( i )*b( i+1, j )-dl( i )*
305 $ b( i+2, j ) ) / d( i )
306 80 CONTINUE
307 IF( j.LT.nrhs ) THEN
308 j = j + 1
309 GO TO 70
310 END IF
311 ELSE
312 DO 100 j = 1, nrhs
313 b( n, j ) = b( n, j ) / d( n )
314 IF( n.GT.1 )
315 $ b( n-1, j ) = ( b( n-1, j )-du( n-1 )*b( n, j ) ) /
316 $ d( n-1 )
317 DO 90 i = n - 2, 1, -1
318 b( i, j ) = ( b( i, j )-du( i )*b( i+1, j )-dl( i )*
319 $ b( i+2, j ) ) / d( i )
320 90 CONTINUE
321 100 CONTINUE
322 END IF
323*
324 RETURN
325*
326* End of DGTSV
327*
328 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgtsv(n, nrhs, dl, d, du, b, ldb, info)
DGTSV computes the solution to system of linear equations A * X = B for GT matrices
Definition dgtsv.f:125