LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sgtrfs.f
Go to the documentation of this file.
1*> \brief \b SGTRFS
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download SGTRFS + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgtrfs.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgtrfs.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgtrfs.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE SGTRFS( TRANS, N, NRHS, DL, D, DU, DLF, DF, DUF, DU2,
22* IPIV, B, LDB, X, LDX, FERR, BERR, WORK, IWORK,
23* INFO )
24*
25* .. Scalar Arguments ..
26* CHARACTER TRANS
27* INTEGER INFO, LDB, LDX, N, NRHS
28* ..
29* .. Array Arguments ..
30* INTEGER IPIV( * ), IWORK( * )
31* REAL B( LDB, * ), BERR( * ), D( * ), DF( * ),
32* $ DL( * ), DLF( * ), DU( * ), DU2( * ), DUF( * ),
33* $ FERR( * ), WORK( * ), X( LDX, * )
34* ..
35*
36*
37*> \par Purpose:
38* =============
39*>
40*> \verbatim
41*>
42*> SGTRFS improves the computed solution to a system of linear
43*> equations when the coefficient matrix is tridiagonal, and provides
44*> error bounds and backward error estimates for the solution.
45*> \endverbatim
46*
47* Arguments:
48* ==========
49*
50*> \param[in] TRANS
51*> \verbatim
52*> TRANS is CHARACTER*1
53*> Specifies the form of the system of equations:
54*> = 'N': A * X = B (No transpose)
55*> = 'T': A**T * X = B (Transpose)
56*> = 'C': A**H * X = B (Conjugate transpose = Transpose)
57*> \endverbatim
58*>
59*> \param[in] N
60*> \verbatim
61*> N is INTEGER
62*> The order of the matrix A. N >= 0.
63*> \endverbatim
64*>
65*> \param[in] NRHS
66*> \verbatim
67*> NRHS is INTEGER
68*> The number of right hand sides, i.e., the number of columns
69*> of the matrix B. NRHS >= 0.
70*> \endverbatim
71*>
72*> \param[in] DL
73*> \verbatim
74*> DL is REAL array, dimension (N-1)
75*> The (n-1) subdiagonal elements of A.
76*> \endverbatim
77*>
78*> \param[in] D
79*> \verbatim
80*> D is REAL array, dimension (N)
81*> The diagonal elements of A.
82*> \endverbatim
83*>
84*> \param[in] DU
85*> \verbatim
86*> DU is REAL array, dimension (N-1)
87*> The (n-1) superdiagonal elements of A.
88*> \endverbatim
89*>
90*> \param[in] DLF
91*> \verbatim
92*> DLF is REAL array, dimension (N-1)
93*> The (n-1) multipliers that define the matrix L from the
94*> LU factorization of A as computed by SGTTRF.
95*> \endverbatim
96*>
97*> \param[in] DF
98*> \verbatim
99*> DF is REAL array, dimension (N)
100*> The n diagonal elements of the upper triangular matrix U from
101*> the LU factorization of A.
102*> \endverbatim
103*>
104*> \param[in] DUF
105*> \verbatim
106*> DUF is REAL array, dimension (N-1)
107*> The (n-1) elements of the first superdiagonal of U.
108*> \endverbatim
109*>
110*> \param[in] DU2
111*> \verbatim
112*> DU2 is REAL array, dimension (N-2)
113*> The (n-2) elements of the second superdiagonal of U.
114*> \endverbatim
115*>
116*> \param[in] IPIV
117*> \verbatim
118*> IPIV is INTEGER array, dimension (N)
119*> The pivot indices; for 1 <= i <= n, row i of the matrix was
120*> interchanged with row IPIV(i). IPIV(i) will always be either
121*> i or i+1; IPIV(i) = i indicates a row interchange was not
122*> required.
123*> \endverbatim
124*>
125*> \param[in] B
126*> \verbatim
127*> B is REAL array, dimension (LDB,NRHS)
128*> The right hand side matrix B.
129*> \endverbatim
130*>
131*> \param[in] LDB
132*> \verbatim
133*> LDB is INTEGER
134*> The leading dimension of the array B. LDB >= max(1,N).
135*> \endverbatim
136*>
137*> \param[in,out] X
138*> \verbatim
139*> X is REAL array, dimension (LDX,NRHS)
140*> On entry, the solution matrix X, as computed by SGTTRS.
141*> On exit, the improved solution matrix X.
142*> \endverbatim
143*>
144*> \param[in] LDX
145*> \verbatim
146*> LDX is INTEGER
147*> The leading dimension of the array X. LDX >= max(1,N).
148*> \endverbatim
149*>
150*> \param[out] FERR
151*> \verbatim
152*> FERR is REAL array, dimension (NRHS)
153*> The estimated forward error bound for each solution vector
154*> X(j) (the j-th column of the solution matrix X).
155*> If XTRUE is the true solution corresponding to X(j), FERR(j)
156*> is an estimated upper bound for the magnitude of the largest
157*> element in (X(j) - XTRUE) divided by the magnitude of the
158*> largest element in X(j). The estimate is as reliable as
159*> the estimate for RCOND, and is almost always a slight
160*> overestimate of the true error.
161*> \endverbatim
162*>
163*> \param[out] BERR
164*> \verbatim
165*> BERR is REAL array, dimension (NRHS)
166*> The componentwise relative backward error of each solution
167*> vector X(j) (i.e., the smallest relative change in
168*> any element of A or B that makes X(j) an exact solution).
169*> \endverbatim
170*>
171*> \param[out] WORK
172*> \verbatim
173*> WORK is REAL array, dimension (3*N)
174*> \endverbatim
175*>
176*> \param[out] IWORK
177*> \verbatim
178*> IWORK is INTEGER array, dimension (N)
179*> \endverbatim
180*>
181*> \param[out] INFO
182*> \verbatim
183*> INFO is INTEGER
184*> = 0: successful exit
185*> < 0: if INFO = -i, the i-th argument had an illegal value
186*> \endverbatim
187*
188*> \par Internal Parameters:
189* =========================
190*>
191*> \verbatim
192*> ITMAX is the maximum number of steps of iterative refinement.
193*> \endverbatim
194*
195* Authors:
196* ========
197*
198*> \author Univ. of Tennessee
199*> \author Univ. of California Berkeley
200*> \author Univ. of Colorado Denver
201*> \author NAG Ltd.
202*
203*> \ingroup gtrfs
204*
205* =====================================================================
206 SUBROUTINE sgtrfs( TRANS, N, NRHS, DL, D, DU, DLF, DF, DUF, DU2,
207 $ IPIV, B, LDB, X, LDX, FERR, BERR, WORK, IWORK,
208 $ INFO )
209*
210* -- LAPACK computational routine --
211* -- LAPACK is a software package provided by Univ. of Tennessee, --
212* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
213*
214* .. Scalar Arguments ..
215 CHARACTER TRANS
216 INTEGER INFO, LDB, LDX, N, NRHS
217* ..
218* .. Array Arguments ..
219 INTEGER IPIV( * ), IWORK( * )
220 REAL B( LDB, * ), BERR( * ), D( * ), DF( * ),
221 $ dl( * ), dlf( * ), du( * ), du2( * ), duf( * ),
222 $ ferr( * ), work( * ), x( ldx, * )
223* ..
224*
225* =====================================================================
226*
227* .. Parameters ..
228 INTEGER ITMAX
229 PARAMETER ( ITMAX = 5 )
230 REAL ZERO, ONE
231 parameter( zero = 0.0e+0, one = 1.0e+0 )
232 REAL TWO
233 parameter( two = 2.0e+0 )
234 REAL THREE
235 parameter( three = 3.0e+0 )
236* ..
237* .. Local Scalars ..
238 LOGICAL NOTRAN
239 CHARACTER TRANSN, TRANST
240 INTEGER COUNT, I, J, KASE, NZ
241 REAL EPS, LSTRES, S, SAFE1, SAFE2, SAFMIN
242* ..
243* .. Local Arrays ..
244 INTEGER ISAVE( 3 )
245* ..
246* .. External Subroutines ..
247 EXTERNAL saxpy, scopy, sgttrs, slacn2, slagtm, xerbla
248* ..
249* .. Intrinsic Functions ..
250 INTRINSIC abs, max
251* ..
252* .. External Functions ..
253 LOGICAL LSAME
254 REAL SLAMCH
255 EXTERNAL lsame, slamch
256* ..
257* .. Executable Statements ..
258*
259* Test the input parameters.
260*
261 info = 0
262 notran = lsame( trans, 'N' )
263 IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
264 $ lsame( trans, 'C' ) ) THEN
265 info = -1
266 ELSE IF( n.LT.0 ) THEN
267 info = -2
268 ELSE IF( nrhs.LT.0 ) THEN
269 info = -3
270 ELSE IF( ldb.LT.max( 1, n ) ) THEN
271 info = -13
272 ELSE IF( ldx.LT.max( 1, n ) ) THEN
273 info = -15
274 END IF
275 IF( info.NE.0 ) THEN
276 CALL xerbla( 'SGTRFS', -info )
277 RETURN
278 END IF
279*
280* Quick return if possible
281*
282 IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
283 DO 10 j = 1, nrhs
284 ferr( j ) = zero
285 berr( j ) = zero
286 10 CONTINUE
287 RETURN
288 END IF
289*
290 IF( notran ) THEN
291 transn = 'N'
292 transt = 'T'
293 ELSE
294 transn = 'T'
295 transt = 'N'
296 END IF
297*
298* NZ = maximum number of nonzero elements in each row of A, plus 1
299*
300 nz = 4
301 eps = slamch( 'Epsilon' )
302 safmin = slamch( 'Safe minimum' )
303 safe1 = nz*safmin
304 safe2 = safe1 / eps
305*
306* Do for each right hand side
307*
308 DO 110 j = 1, nrhs
309*
310 count = 1
311 lstres = three
312 20 CONTINUE
313*
314* Loop until stopping criterion is satisfied.
315*
316* Compute residual R = B - op(A) * X,
317* where op(A) = A, A**T, or A**H, depending on TRANS.
318*
319 CALL scopy( n, b( 1, j ), 1, work( n+1 ), 1 )
320 CALL slagtm( trans, n, 1, -one, dl, d, du, x( 1, j ), ldx, one,
321 $ work( n+1 ), n )
322*
323* Compute abs(op(A))*abs(x) + abs(b) for use in the backward
324* error bound.
325*
326 IF( notran ) THEN
327 IF( n.EQ.1 ) THEN
328 work( 1 ) = abs( b( 1, j ) ) + abs( d( 1 )*x( 1, j ) )
329 ELSE
330 work( 1 ) = abs( b( 1, j ) ) + abs( d( 1 )*x( 1, j ) ) +
331 $ abs( du( 1 )*x( 2, j ) )
332 DO 30 i = 2, n - 1
333 work( i ) = abs( b( i, j ) ) +
334 $ abs( dl( i-1 )*x( i-1, j ) ) +
335 $ abs( d( i )*x( i, j ) ) +
336 $ abs( du( i )*x( i+1, j ) )
337 30 CONTINUE
338 work( n ) = abs( b( n, j ) ) +
339 $ abs( dl( n-1 )*x( n-1, j ) ) +
340 $ abs( d( n )*x( n, j ) )
341 END IF
342 ELSE
343 IF( n.EQ.1 ) THEN
344 work( 1 ) = abs( b( 1, j ) ) + abs( d( 1 )*x( 1, j ) )
345 ELSE
346 work( 1 ) = abs( b( 1, j ) ) + abs( d( 1 )*x( 1, j ) ) +
347 $ abs( dl( 1 )*x( 2, j ) )
348 DO 40 i = 2, n - 1
349 work( i ) = abs( b( i, j ) ) +
350 $ abs( du( i-1 )*x( i-1, j ) ) +
351 $ abs( d( i )*x( i, j ) ) +
352 $ abs( dl( i )*x( i+1, j ) )
353 40 CONTINUE
354 work( n ) = abs( b( n, j ) ) +
355 $ abs( du( n-1 )*x( n-1, j ) ) +
356 $ abs( d( n )*x( n, j ) )
357 END IF
358 END IF
359*
360* Compute componentwise relative backward error from formula
361*
362* max(i) ( abs(R(i)) / ( abs(op(A))*abs(X) + abs(B) )(i) )
363*
364* where abs(Z) is the componentwise absolute value of the matrix
365* or vector Z. If the i-th component of the denominator is less
366* than SAFE2, then SAFE1 is added to the i-th components of the
367* numerator and denominator before dividing.
368*
369 s = zero
370 DO 50 i = 1, n
371 IF( work( i ).GT.safe2 ) THEN
372 s = max( s, abs( work( n+i ) ) / work( i ) )
373 ELSE
374 s = max( s, ( abs( work( n+i ) )+safe1 ) /
375 $ ( work( i )+safe1 ) )
376 END IF
377 50 CONTINUE
378 berr( j ) = s
379*
380* Test stopping criterion. Continue iterating if
381* 1) The residual BERR(J) is larger than machine epsilon, and
382* 2) BERR(J) decreased by at least a factor of 2 during the
383* last iteration, and
384* 3) At most ITMAX iterations tried.
385*
386 IF( berr( j ).GT.eps .AND. two*berr( j ).LE.lstres .AND.
387 $ count.LE.itmax ) THEN
388*
389* Update solution and try again.
390*
391 CALL sgttrs( trans, n, 1, dlf, df, duf, du2, ipiv,
392 $ work( n+1 ), n, info )
393 CALL saxpy( n, one, work( n+1 ), 1, x( 1, j ), 1 )
394 lstres = berr( j )
395 count = count + 1
396 GO TO 20
397 END IF
398*
399* Bound error from formula
400*
401* norm(X - XTRUE) / norm(X) .le. FERR =
402* norm( abs(inv(op(A)))*
403* ( abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) ))) / norm(X)
404*
405* where
406* norm(Z) is the magnitude of the largest component of Z
407* inv(op(A)) is the inverse of op(A)
408* abs(Z) is the componentwise absolute value of the matrix or
409* vector Z
410* NZ is the maximum number of nonzeros in any row of A, plus 1
411* EPS is machine epsilon
412*
413* The i-th component of abs(R)+NZ*EPS*(abs(op(A))*abs(X)+abs(B))
414* is incremented by SAFE1 if the i-th component of
415* abs(op(A))*abs(X) + abs(B) is less than SAFE2.
416*
417* Use SLACN2 to estimate the infinity-norm of the matrix
418* inv(op(A)) * diag(W),
419* where W = abs(R) + NZ*EPS*( abs(op(A))*abs(X)+abs(B) )))
420*
421 DO 60 i = 1, n
422 IF( work( i ).GT.safe2 ) THEN
423 work( i ) = abs( work( n+i ) ) + nz*eps*work( i )
424 ELSE
425 work( i ) = abs( work( n+i ) ) + nz*eps*work( i ) + safe1
426 END IF
427 60 CONTINUE
428*
429 kase = 0
430 70 CONTINUE
431 CALL slacn2( n, work( 2*n+1 ), work( n+1 ), iwork, ferr( j ),
432 $ kase, isave )
433 IF( kase.NE.0 ) THEN
434 IF( kase.EQ.1 ) THEN
435*
436* Multiply by diag(W)*inv(op(A)**T).
437*
438 CALL sgttrs( transt, n, 1, dlf, df, duf, du2, ipiv,
439 $ work( n+1 ), n, info )
440 DO 80 i = 1, n
441 work( n+i ) = work( i )*work( n+i )
442 80 CONTINUE
443 ELSE
444*
445* Multiply by inv(op(A))*diag(W).
446*
447 DO 90 i = 1, n
448 work( n+i ) = work( i )*work( n+i )
449 90 CONTINUE
450 CALL sgttrs( transn, n, 1, dlf, df, duf, du2, ipiv,
451 $ work( n+1 ), n, info )
452 END IF
453 GO TO 70
454 END IF
455*
456* Normalize error.
457*
458 lstres = zero
459 DO 100 i = 1, n
460 lstres = max( lstres, abs( x( i, j ) ) )
461 100 CONTINUE
462 IF( lstres.NE.zero )
463 $ ferr( j ) = ferr( j ) / lstres
464*
465 110 CONTINUE
466*
467 RETURN
468*
469* End of SGTRFS
470*
471 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine saxpy(n, sa, sx, incx, sy, incy)
SAXPY
Definition saxpy.f:89
subroutine scopy(n, sx, incx, sy, incy)
SCOPY
Definition scopy.f:82
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:209
subroutine sgttrs(trans, n, nrhs, dl, d, du, du2, ipiv, b, ldb, info)
SGTTRS
Definition sgttrs.f:138
subroutine slacn2(n, v, x, isgn, est, kase, isave)
SLACN2 estimates the 1-norm of a square matrix, using reverse communication for evaluating matrix-vec...
Definition slacn2.f:136
subroutine slagtm(trans, n, nrhs, alpha, dl, d, du, x, ldx, beta, b, ldb)
SLAGTM performs a matrix-matrix product of the form C = αAB+βC, where A is a tridiagonal matrix,...
Definition slagtm.f:145