LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sgtsvx.f
Go to the documentation of this file.
1*> \brief <b> SGTSVX 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 SGTSVX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgtsvx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgtsvx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgtsvx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SGTSVX( FACT, TRANS, N, NRHS, DL, D, DU, DLF, DF, DUF,
20* DU2, IPIV, B, LDB, X, LDX, RCOND, FERR, BERR,
21* WORK, IWORK, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER FACT, TRANS
25* INTEGER INFO, LDB, LDX, N, NRHS
26* REAL RCOND
27* ..
28* .. Array Arguments ..
29* INTEGER IPIV( * ), IWORK( * )
30* REAL B( LDB, * ), BERR( * ), D( * ), DF( * ),
31* $ DL( * ), DLF( * ), DU( * ), DU2( * ), DUF( * ),
32* $ FERR( * ), WORK( * ), X( LDX, * )
33* ..
34*
35*
36*> \par Purpose:
37* =============
38*>
39*> \verbatim
40*>
41*> SGTSVX uses the LU factorization to compute the solution to a real
42*> system of linear equations A * X = B or A**T * X = B,
43*> where A is a tridiagonal matrix of order N and X and B are N-by-NRHS
44*> matrices.
45*>
46*> Error bounds on the solution and a condition estimate are also
47*> provided.
48*> \endverbatim
49*
50*> \par Description:
51* =================
52*>
53*> \verbatim
54*>
55*> The following steps are performed:
56*>
57*> 1. If FACT = 'N', the LU decomposition is used to factor the matrix A
58*> as A = L * U, where L is a product of permutation and unit lower
59*> bidiagonal matrices and U is upper triangular with nonzeros in
60*> only the main diagonal and first two superdiagonals.
61*>
62*> 2. If some U(i,i)=0, so that U is exactly singular, then the routine
63*> returns with INFO = i. Otherwise, the factored form of A is used
64*> to estimate the condition number of the matrix A. If the
65*> reciprocal of the condition number is less than machine precision,
66*> INFO = N+1 is returned as a warning, but the routine still goes on
67*> to solve for X and compute error bounds as described below.
68*>
69*> 3. The system of equations is solved for X using the factored form
70*> of A.
71*>
72*> 4. Iterative refinement is applied to improve the computed solution
73*> matrix and calculate error bounds and backward error estimates
74*> for it.
75*> \endverbatim
76*
77* Arguments:
78* ==========
79*
80*> \param[in] FACT
81*> \verbatim
82*> FACT is CHARACTER*1
83*> Specifies whether or not the factored form of A has been
84*> supplied on entry.
85*> = 'F': DLF, DF, DUF, DU2, and IPIV contain the factored
86*> form of A; DL, D, DU, DLF, DF, DUF, DU2 and IPIV
87*> will not be modified.
88*> = 'N': The matrix will be copied to DLF, DF, and DUF
89*> and factored.
90*> \endverbatim
91*>
92*> \param[in] TRANS
93*> \verbatim
94*> TRANS is CHARACTER*1
95*> Specifies the form of the system of equations:
96*> = 'N': A * X = B (No transpose)
97*> = 'T': A**T * X = B (Transpose)
98*> = 'C': A**H * X = B (Conjugate transpose = Transpose)
99*> \endverbatim
100*>
101*> \param[in] N
102*> \verbatim
103*> N is INTEGER
104*> The order of the matrix A. N >= 0.
105*> \endverbatim
106*>
107*> \param[in] NRHS
108*> \verbatim
109*> NRHS is INTEGER
110*> The number of right hand sides, i.e., the number of columns
111*> of the matrix B. NRHS >= 0.
112*> \endverbatim
113*>
114*> \param[in] DL
115*> \verbatim
116*> DL is REAL array, dimension (N-1)
117*> The (n-1) subdiagonal elements of A.
118*> \endverbatim
119*>
120*> \param[in] D
121*> \verbatim
122*> D is REAL array, dimension (N)
123*> The n diagonal elements of A.
124*> \endverbatim
125*>
126*> \param[in] DU
127*> \verbatim
128*> DU is REAL array, dimension (N-1)
129*> The (n-1) superdiagonal elements of A.
130*> \endverbatim
131*>
132*> \param[in,out] DLF
133*> \verbatim
134*> DLF is REAL array, dimension (N-1)
135*> If FACT = 'F', then DLF is an input argument and on entry
136*> contains the (n-1) multipliers that define the matrix L from
137*> the LU factorization of A as computed by SGTTRF.
138*>
139*> If FACT = 'N', then DLF is an output argument and on exit
140*> contains the (n-1) multipliers that define the matrix L from
141*> the LU factorization of A.
142*> \endverbatim
143*>
144*> \param[in,out] DF
145*> \verbatim
146*> DF is REAL array, dimension (N)
147*> If FACT = 'F', then DF is an input argument and on entry
148*> contains the n diagonal elements of the upper triangular
149*> matrix U from the LU factorization of A.
150*>
151*> If FACT = 'N', then DF is an output argument and on exit
152*> contains the n diagonal elements of the upper triangular
153*> matrix U from the LU factorization of A.
154*> \endverbatim
155*>
156*> \param[in,out] DUF
157*> \verbatim
158*> DUF is REAL array, dimension (N-1)
159*> If FACT = 'F', then DUF is an input argument and on entry
160*> contains the (n-1) elements of the first superdiagonal of U.
161*>
162*> If FACT = 'N', then DUF is an output argument and on exit
163*> contains the (n-1) elements of the first superdiagonal of U.
164*> \endverbatim
165*>
166*> \param[in,out] DU2
167*> \verbatim
168*> DU2 is REAL array, dimension (N-2)
169*> If FACT = 'F', then DU2 is an input argument and on entry
170*> contains the (n-2) elements of the second superdiagonal of
171*> U.
172*>
173*> If FACT = 'N', then DU2 is an output argument and on exit
174*> contains the (n-2) elements of the second superdiagonal of
175*> U.
176*> \endverbatim
177*>
178*> \param[in,out] IPIV
179*> \verbatim
180*> IPIV is INTEGER array, dimension (N)
181*> If FACT = 'F', then IPIV is an input argument and on entry
182*> contains the pivot indices from the LU factorization of A as
183*> computed by SGTTRF.
184*>
185*> If FACT = 'N', then IPIV is an output argument and on exit
186*> contains the pivot indices from the LU factorization of A;
187*> row i of the matrix was interchanged with row IPIV(i).
188*> IPIV(i) will always be either i or i+1; IPIV(i) = i indicates
189*> a row interchange was not required.
190*> \endverbatim
191*>
192*> \param[in] B
193*> \verbatim
194*> B is REAL array, dimension (LDB,NRHS)
195*> The N-by-NRHS right hand side matrix B.
196*> \endverbatim
197*>
198*> \param[in] LDB
199*> \verbatim
200*> LDB is INTEGER
201*> The leading dimension of the array B. LDB >= max(1,N).
202*> \endverbatim
203*>
204*> \param[out] X
205*> \verbatim
206*> X is REAL array, dimension (LDX,NRHS)
207*> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X.
208*> \endverbatim
209*>
210*> \param[in] LDX
211*> \verbatim
212*> LDX is INTEGER
213*> The leading dimension of the array X. LDX >= max(1,N).
214*> \endverbatim
215*>
216*> \param[out] RCOND
217*> \verbatim
218*> RCOND is REAL
219*> The estimate of the reciprocal condition number of the matrix
220*> A. If RCOND is less than the machine precision (in
221*> particular, if RCOND = 0), the matrix is singular to working
222*> precision. This condition is indicated by a return code of
223*> INFO > 0.
224*> \endverbatim
225*>
226*> \param[out] FERR
227*> \verbatim
228*> FERR is REAL array, dimension (NRHS)
229*> The estimated forward error bound for each solution vector
230*> X(j) (the j-th column of the solution matrix X).
231*> If XTRUE is the true solution corresponding to X(j), FERR(j)
232*> is an estimated upper bound for the magnitude of the largest
233*> element in (X(j) - XTRUE) divided by the magnitude of the
234*> largest element in X(j). The estimate is as reliable as
235*> the estimate for RCOND, and is almost always a slight
236*> overestimate of the true error.
237*> \endverbatim
238*>
239*> \param[out] BERR
240*> \verbatim
241*> BERR is REAL array, dimension (NRHS)
242*> The componentwise relative backward error of each solution
243*> vector X(j) (i.e., the smallest relative change in
244*> any element of A or B that makes X(j) an exact solution).
245*> \endverbatim
246*>
247*> \param[out] WORK
248*> \verbatim
249*> WORK is REAL array, dimension (3*N)
250*> \endverbatim
251*>
252*> \param[out] IWORK
253*> \verbatim
254*> IWORK is INTEGER array, dimension (N)
255*> \endverbatim
256*>
257*> \param[out] INFO
258*> \verbatim
259*> INFO is INTEGER
260*> = 0: successful exit
261*> < 0: if INFO = -i, the i-th argument had an illegal value
262*> > 0: if INFO = i, and i is
263*> <= N: U(i,i) is exactly zero. The factorization
264*> has not been completed unless i = N, but the
265*> factor U is exactly singular, so the solution
266*> and error bounds could not be computed.
267*> RCOND = 0 is returned.
268*> = N+1: U is nonsingular, but RCOND is less than machine
269*> precision, meaning that the matrix is singular
270*> to working precision. Nevertheless, the
271*> solution and error bounds are computed because
272*> there are a number of situations where the
273*> computed solution can be more accurate than the
274*> value of RCOND would suggest.
275*> \endverbatim
276*
277* Authors:
278* ========
279*
280*> \author Univ. of Tennessee
281*> \author Univ. of California Berkeley
282*> \author Univ. of Colorado Denver
283*> \author NAG Ltd.
284*
285*> \ingroup gtsvx
286*
287* =====================================================================
288 SUBROUTINE sgtsvx( FACT, TRANS, N, NRHS, DL, D, DU, DLF, DF,
289 $ DUF,
290 $ DU2, IPIV, B, LDB, X, LDX, RCOND, FERR, BERR,
291 $ WORK, IWORK, INFO )
292*
293* -- LAPACK driver routine --
294* -- LAPACK is a software package provided by Univ. of Tennessee, --
295* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
296*
297* .. Scalar Arguments ..
298 CHARACTER FACT, TRANS
299 INTEGER INFO, LDB, LDX, N, NRHS
300 REAL RCOND
301* ..
302* .. Array Arguments ..
303 INTEGER IPIV( * ), IWORK( * )
304 REAL B( LDB, * ), BERR( * ), D( * ), DF( * ),
305 $ DL( * ), DLF( * ), DU( * ), DU2( * ), DUF( * ),
306 $ ferr( * ), work( * ), x( ldx, * )
307* ..
308*
309* =====================================================================
310*
311* .. Parameters ..
312 REAL ZERO
313 PARAMETER ( ZERO = 0.0e+0 )
314* ..
315* .. Local Scalars ..
316 LOGICAL NOFACT, NOTRAN
317 CHARACTER NORM
318 REAL ANORM
319* ..
320* .. External Functions ..
321 LOGICAL LSAME
322 REAL SLAMCH, SLANGT
323 EXTERNAL LSAME, SLAMCH, SLANGT
324* ..
325* .. External Subroutines ..
326 EXTERNAL scopy, sgtcon, sgtrfs, sgttrf, sgttrs,
327 $ slacpy,
328 $ xerbla
329* ..
330* .. Intrinsic Functions ..
331 INTRINSIC max
332* ..
333* .. Executable Statements ..
334*
335 info = 0
336 nofact = lsame( fact, 'N' )
337 notran = lsame( trans, 'N' )
338 IF( .NOT.nofact .AND. .NOT.lsame( fact, 'F' ) ) THEN
339 info = -1
340 ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
341 $ lsame( trans, 'C' ) ) THEN
342 info = -2
343 ELSE IF( n.LT.0 ) THEN
344 info = -3
345 ELSE IF( nrhs.LT.0 ) THEN
346 info = -4
347 ELSE IF( ldb.LT.max( 1, n ) ) THEN
348 info = -14
349 ELSE IF( ldx.LT.max( 1, n ) ) THEN
350 info = -16
351 END IF
352 IF( info.NE.0 ) THEN
353 CALL xerbla( 'SGTSVX', -info )
354 RETURN
355 END IF
356*
357 IF( nofact ) THEN
358*
359* Compute the LU factorization of A.
360*
361 CALL scopy( n, d, 1, df, 1 )
362 IF( n.GT.1 ) THEN
363 CALL scopy( n-1, dl, 1, dlf, 1 )
364 CALL scopy( n-1, du, 1, duf, 1 )
365 END IF
366 CALL sgttrf( n, dlf, df, duf, du2, ipiv, info )
367*
368* Return if INFO is non-zero.
369*
370 IF( info.GT.0 )THEN
371 rcond = zero
372 RETURN
373 END IF
374 END IF
375*
376* Compute the norm of the matrix A.
377*
378 IF( notran ) THEN
379 norm = '1'
380 ELSE
381 norm = 'I'
382 END IF
383 anorm = slangt( norm, n, dl, d, du )
384*
385* Compute the reciprocal of the condition number of A.
386*
387 CALL sgtcon( norm, n, dlf, df, duf, du2, ipiv, anorm, rcond,
388 $ work,
389 $ iwork, info )
390*
391* Compute the solution vectors X.
392*
393 CALL slacpy( 'Full', n, nrhs, b, ldb, x, ldx )
394 CALL sgttrs( trans, n, nrhs, dlf, df, duf, du2, ipiv, x, ldx,
395 $ info )
396*
397* Use iterative refinement to improve the computed solutions and
398* compute error bounds and backward error estimates for them.
399*
400 CALL sgtrfs( trans, n, nrhs, dl, d, du, dlf, df, duf, du2,
401 $ ipiv,
402 $ b, ldb, x, ldx, ferr, berr, work, iwork, info )
403*
404* Set INFO = N+1 if the matrix is singular to working precision.
405*
406 IF( rcond.LT.slamch( 'Epsilon' ) )
407 $ info = n + 1
408*
409 RETURN
410*
411* End of SGTSVX
412*
413 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
subroutine sgtcon(norm, n, dl, d, du, du2, ipiv, anorm, rcond, work, iwork, info)
SGTCON
Definition sgtcon.f:144
subroutine sgtrfs(trans, n, nrhs, dl, d, du, dlf, df, duf, du2, ipiv, b, ldb, x, ldx, ferr, berr, work, iwork, info)
SGTRFS
Definition sgtrfs.f:208
subroutine sgtsvx(fact, trans, n, nrhs, dl, d, du, dlf, df, duf, du2, ipiv, b, ldb, x, ldx, rcond, ferr, berr, work, iwork, info)
SGTSVX computes the solution to system of linear equations A * X = B for GT matrices
Definition sgtsvx.f:292
subroutine sgttrf(n, dl, d, du, du2, ipiv, info)
SGTTRF
Definition sgttrf.f:122
subroutine sgttrs(trans, n, nrhs, dl, d, du, du2, ipiv, b, ldb, info)
SGTTRS
Definition sgttrs.f:137
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:101