1SGBSVXX(1) LAPACK driver routine (version 3.2) SGBSVXX(1)
2
3
4
6 SGBSVXX - SGBSVXX use the LU factorization to compute the solution to a
7 real system of linear equations A * X = B, where A is an N-by-N matrix
8 and X and B are N-by-NRHS matrices
9
11 SUBROUTINE SGBSVXX( FACT, TRANS, N, KL, KU, NRHS, AB, LDAB, AFB, LDAFB,
12 IPIV, EQUED, R, C, B, LDB, X, LDX, RCOND, RPVGRW,
13 BERR, N_ERR_BNDS, ERR_BNDS_NORM, ERR_BNDS_COMP,
14 NPARAMS, PARAMS, WORK, IWORK, INFO )
15
16 IMPLICIT NONE
17
18 CHARACTER EQUED, FACT, TRANS
19
20 INTEGER INFO, LDAB, LDAFB, LDB, LDX, N, NRHS, NPARAMS,
21 N_ERR_BNDS
22
23 REAL RCOND, RPVGRW
24
25 INTEGER IPIV( * ), IWORK( * )
26
27 REAL AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ), X( LDX
28 , * ),WORK( * )
29
30 REAL R( * ), C( * ), PARAMS( * ), BERR( * ),
31 ERR_BNDS_NORM( NRHS, * ), ERR_BNDS_COMP( NRHS, * )
32
34 SGBSVXX uses the LU factorization to compute the solution to a
35 real system of linear equations A * X = B, where A is an
36 N-by-N matrix and X and B are N-by-NRHS matrices.
37 If requested, both normwise and maximum componentwise error bounds
38 are returned. SGBSVXX will return a solution with a tiny
39 guaranteed error (O(eps) where eps is the working machine
40 precision) unless the matrix is very ill-conditioned, in which
41 case a warning is returned. Relevant condition numbers also are
42 calculated and returned.
43 SGBSVXX accepts user-provided factorizations and equilibration
44 factors; see the definitions of the FACT and EQUED options.
45 Solving with refinement and using a factorization from a previous
46 SGBSVXX call will also produce a solution with either O(eps)
47 errors or warnings, but we cannot make that claim for general
48 user-provided factorizations and equilibration factors if they
49 differ from what SGBSVXX would itself produce.
50
52 The following steps are performed:
53 1. If FACT = 'E', real scaling factors are computed to equilibrate
54 the system:
55 TRANS = 'N': diag(R)*A*diag(C) *inv(diag(C))*X = diag(R)*B
56 TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
57 TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B
58 Whether or not the system will be equilibrated depends on the
59 scaling of the matrix A, but if equilibration is used, A is
60 overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N')
61 or diag(C)*B (if TRANS = 'T' or 'C').
62 2. If FACT = 'N' or 'E', the LU decomposition is used to factor
63 the matrix A (after equilibration if FACT = 'E') as
64 A = P * L * U,
65 where P is a permutation matrix, L is a unit lower triangular
66 matrix, and U is upper triangular.
67 3. If some U(i,i)=0, so that U is exactly singular, then the
68 routine returns with INFO = i. Otherwise, the factored form of A
69 is used to estimate the condition number of the matrix A (see
70 argument RCOND). If the reciprocal of the condition number is less
71 than machine precision, the routine still goes on to solve for X
72 and compute error bounds as described below.
73 4. The system of equations is solved for X using the factored form
74 of A.
75 5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
76 the routine will use iterative refinement to try to get a small
77 error and error bounds. Refinement calculates the residual to at
78 least twice the working precision.
79 6. If equilibration was used, the matrix X is premultiplied by
80 diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so
81 that it solves the original system before equilibration.
82
84 Some optional parameters are bundled in the PARAMS array. These set‐
85 tings determine how refinement is performed, but often the defaults are
86 acceptable. If the defaults are acceptable, users can pass NPARAMS = 0
87 which prevents the source code from accessing the PARAMS argument.
88
89 FACT (input) CHARACTER*1
90 Specifies whether or not the factored form of the matrix A is
91 supplied on entry, and if not, whether the matrix A should be
92 equilibrated before it is factored. = 'F': On entry, AF and
93 IPIV contain the factored form of A. If EQUED is not 'N', the
94 matrix A has been equilibrated with scaling factors given by R
95 and C. A, AF, and IPIV are not modified. = 'N': The matrix A
96 will be copied to AF and factored.
97 = 'E': The matrix A will be equilibrated if necessary, then
98 copied to AF and factored.
99
100 TRANS (input) CHARACTER*1
101 Specifies the form of the system of equations:
102 = 'N': A * X = B (No transpose)
103 = 'T': A**T * X = B (Transpose)
104 = 'C': A**H * X = B (Conjugate Transpose = Transpose)
105
106 N (input) INTEGER
107 The number of linear equations, i.e., the order of the matrix
108 A. N >= 0.
109
110 KL (input) INTEGER
111 The number of subdiagonals within the band of A. KL >= 0.
112
113 KU (input) INTEGER
114 The number of superdiagonals within the band of A. KU >= 0.
115
116 NRHS (input) INTEGER
117 The number of right hand sides, i.e., the number of columns of
118 the matrices B and X. NRHS >= 0.
119
120 AB (input/output) REAL array, dimension (LDAB,N)
121 On entry, the matrix A in band storage, in rows 1 to KL+KU+1.
122 The j-th column of A is stored in the j-th column of the array
123 AB as follows: AB(KU+1+i-j,j) = A(i,j) for max(1,j-
124 KU)<=i<=min(N,j+kl) If FACT = 'F' and EQUED is not 'N', then AB
125 must have been equilibrated by the scaling factors in R and/or
126 C. AB is not modified if FACT = 'F' or 'N', or if FACT = 'E'
127 and EQUED = 'N' on exit. On exit, if EQUED .ne. 'N', A is
128 scaled as follows: EQUED = 'R': A := diag(R) * A
129 EQUED = 'C': A := A * diag(C)
130 EQUED = 'B': A := diag(R) * A * diag(C).
131
132 LDAB (input) INTEGER
133 The leading dimension of the array AB. LDAB >= KL+KU+1.
134
135 AFB (input or output) REAL array, dimension (LDAFB,N)
136 If FACT = 'F', then AFB is an input argument and on entry con‐
137 tains details of the LU factorization of the band matrix A, as
138 computed by SGBTRF. U is stored as an upper triangular band
139 matrix with KL+KU superdiagonals in rows 1 to KL+KU+1, and the
140 multipliers used during the factorization are stored in rows
141 KL+KU+2 to 2*KL+KU+1. If EQUED .ne. 'N', then AFB is the fac‐
142 tored form of the equilibrated matrix A. If FACT = 'N', then
143 AF is an output argument and on exit returns the factors L and
144 U from the factorization A = P*L*U of the original matrix A.
145 If FACT = 'E', then AF is an output argument and on exit
146 returns the factors L and U from the factorization A = P*L*U of
147 the equilibrated matrix A (see the description of A for the
148 form of the equilibrated matrix).
149
150 LDAFB (input) INTEGER
151 The leading dimension of the array AFB. LDAFB >= 2*KL+KU+1.
152
153 IPIV (input or output) INTEGER array, dimension (N)
154 If FACT = 'F', then IPIV is an input argument and on entry con‐
155 tains the pivot indices from the factorization A = P*L*U as
156 computed by SGETRF; row i of the matrix was interchanged with
157 row IPIV(i). If FACT = 'N', then IPIV is an output argument
158 and on exit contains the pivot indices from the factorization A
159 = P*L*U of the original matrix A. If FACT = 'E', then IPIV is
160 an output argument and on exit contains the pivot indices from
161 the factorization A = P*L*U of the equilibrated matrix A.
162
163 EQUED (input or output) CHARACTER*1
164 Specifies the form of equilibration that was done. = 'N': No
165 equilibration (always true if FACT = 'N').
166 = 'R': Row equilibration, i.e., A has been premultiplied by
167 diag(R). = 'C': Column equilibration, i.e., A has been post‐
168 multiplied by diag(C). = 'B': Both row and column equilibra‐
169 tion, i.e., A has been replaced by diag(R) * A * diag(C).
170 EQUED is an input argument if FACT = 'F'; otherwise, it is an
171 output argument.
172
173 R (input or output) REAL array, dimension (N)
174 The row scale factors for A. If EQUED = 'R' or 'B', A is mul‐
175 tiplied on the left by diag(R); if EQUED = 'N' or 'C', R is not
176 accessed. R is an input argument if FACT = 'F'; otherwise, R
177 is an output argument. If FACT = 'F' and EQUED = 'R' or 'B',
178 each element of R must be positive. If R is output, each ele‐
179 ment of R is a power of the radix. If R is input, each element
180 of R should be a power of the radix to ensure a reliable solu‐
181 tion and error estimates. Scaling by powers of the radix does
182 not cause rounding errors unless the result underflows or over‐
183 flows. Rounding errors during scaling lead to refining with a
184 matrix that is not equivalent to the input matrix, producing
185 error estimates that may not be reliable.
186
187 C (input or output) REAL array, dimension (N)
188 The column scale factors for A. If EQUED = 'C' or 'B', A is
189 multiplied on the right by diag(C); if EQUED = 'N' or 'R', C is
190 not accessed. C is an input argument if FACT = 'F'; otherwise,
191 C is an output argument. If FACT = 'F' and EQUED = 'C' or 'B',
192 each element of C must be positive. If C is output, each ele‐
193 ment of C is a power of the radix. If C is input, each element
194 of C should be a power of the radix to ensure a reliable solu‐
195 tion and error estimates. Scaling by powers of the radix does
196 not cause rounding errors unless the result underflows or over‐
197 flows. Rounding errors during scaling lead to refining with a
198 matrix that is not equivalent to the input matrix, producing
199 error estimates that may not be reliable.
200
201 B (input/output) REAL array, dimension (LDB,NRHS)
202 On entry, the N-by-NRHS right hand side matrix B. On exit, if
203 EQUED = 'N', B is not modified; if TRANS = 'N' and EQUED = 'R'
204 or 'B', B is overwritten by diag(R)*B; if TRANS = 'T' or 'C'
205 and EQUED = 'C' or 'B', B is overwritten by diag(C)*B.
206
207 LDB (input) INTEGER
208 The leading dimension of the array B. LDB >= max(1,N).
209
210 X (output) REAL array, dimension (LDX,NRHS)
211 If INFO = 0, the N-by-NRHS solution matrix X to the original
212 system of equations. Note that A and B are modified on exit if
213 EQUED .ne. 'N', and the solution to the equilibrated system is
214 inv(diag(C))*X if TRANS = 'N' and EQUED = 'C' or 'B', or
215 inv(diag(R))*X if TRANS = 'T' or 'C' and EQUED = 'R' or 'B'.
216
217 LDX (input) INTEGER
218 The leading dimension of the array X. LDX >= max(1,N).
219
220 RCOND (output) REAL
221 Reciprocal scaled condition number. This is an estimate of the
222 reciprocal Skeel condition number of the matrix A after equili‐
223 bration (if done). If this is less than the machine precision
224 (in particular, if it is zero), the matrix is singular to work‐
225 ing precision. Note that the error may still be small even if
226 this number is very small and the matrix appears ill- condi‐
227 tioned.
228
229 RPVGRW (output) REAL
230 Reciprocal pivot growth. On exit, this contains the reciprocal
231 pivot growth factor norm(A)/norm(U). The "max absolute element"
232 norm is used. If this is much less than 1, then the stability
233 of the LU factorization of the (equilibrated) matrix A could be
234 poor. This also means that the solution X, estimated condition
235 numbers, and error bounds could be unreliable. If factorization
236 fails with 0<INFO<=N, then this contains the reciprocal pivot
237 growth factor for the leading INFO columns of A. In SGESVX,
238 this quantity is returned in WORK(1).
239
240 BERR (output) REAL array, dimension (NRHS)
241 Componentwise relative backward error. This is the component‐
242 wise relative backward error of each solution vector X(j)
243 (i.e., the smallest relative change in any element of A or B
244 that makes X(j) an exact solution). N_ERR_BNDS (input) INTEGER
245 Number of error bounds to return for each right hand side and
246 each type (normwise or componentwise). See ERR_BNDS_NORM and
247 ERR_BNDS_COMP below.
248
249 ERR_BNDS_NORM (output) REAL array, dimension (NRHS, N_ERR_BNDS)
250 For each right-hand side, this array contains informa‐
251 tion about various error bounds and condition numbers
252 corresponding to the normwise relative error, which is
253 defined as follows: Normwise relative error in the ith
254 solution vector: max_j (abs(XTRUE(j,i) - X(j,i)))
255 ------------------------------ max_j abs(X(j,i)) The
256 array is indexed by the type of error information as
257 described below. There currently are up to three pieces
258 of information returned. The first index in
259 ERR_BNDS_NORM(i,:) corresponds to the ith right-hand
260 side. The second index in ERR_BNDS_NORM(:,err) contains
261 the following three fields: err = 1 "Trust/don't trust"
262 boolean. Trust the answer if the reciprocal condition
263 number is less than the threshold sqrt(n) *
264 slamch('Epsilon'). err = 2 "Guaranteed" error bound:
265 The estimated forward error, almost certainly within a
266 factor of 10 of the true error so long as the next entry
267 is greater than the threshold sqrt(n) *
268 slamch('Epsilon'). This error bound should only be
269 trusted if the previous boolean is true. err = 3
270 Reciprocal condition number: Estimated normwise recipro‐
271 cal condition number. Compared with the threshold
272 sqrt(n) * slamch('Epsilon') to determine if the error
273 estimate is "guaranteed". These reciprocal condition
274 numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf)) for
275 some appropriately scaled matrix Z. Let Z = S*A, where
276 S scales each row by a power of the radix so all abso‐
277 lute row sums of Z are approximately 1. See Lapack
278 Working Note 165 for further details and extra cautions.
279
280 ERR_BNDS_COMP (output) REAL array, dimension (NRHS, N_ERR_BNDS)
281 For each right-hand side, this array contains informa‐
282 tion about various error bounds and condition numbers
283 corresponding to the componentwise relative error, which
284 is defined as follows: Componentwise relative error in
285 the ith solution vector: abs(XTRUE(j,i) - X(j,i)) max_j
286 ---------------------- abs(X(j,i)) The array is indexed
287 by the right-hand side i (on which the componentwise
288 relative error depends), and the type of error informa‐
289 tion as described below. There currently are up to three
290 pieces of information returned for each right-hand side.
291 If componentwise accuracy is not requested (PARAMS(3) =
292 0.0), then ERR_BNDS_COMP is not accessed. If N_ERR_BNDS
293 .LT. 3, then at most the first (:,N_ERR_BNDS) entries
294 are returned. The first index in ERR_BNDS_COMP(i,:)
295 corresponds to the ith right-hand side. The second
296 index in ERR_BNDS_COMP(:,err) contains the following
297 three fields: err = 1 "Trust/don't trust" boolean. Trust
298 the answer if the reciprocal condition number is less
299 than the threshold sqrt(n) * slamch('Epsilon'). err = 2
300 "Guaranteed" error bound: The estimated forward error,
301 almost certainly within a factor of 10 of the true error
302 so long as the next entry is greater than the threshold
303 sqrt(n) * slamch('Epsilon'). This error bound should
304 only be trusted if the previous boolean is true. err =
305 3 Reciprocal condition number: Estimated componentwise
306 reciprocal condition number. Compared with the thresh‐
307 old sqrt(n) * slamch('Epsilon') to determine if the
308 error estimate is "guaranteed". These reciprocal condi‐
309 tion numbers are 1 / (norm(Z^{-1},inf) * norm(Z,inf))
310 for some appropriately scaled matrix Z. Let Z =
311 S*(A*diag(x)), where x is the solution for the current
312 right-hand side and S scales each row of A*diag(x) by a
313 power of the radix so all absolute row sums of Z are
314 approximately 1. See Lapack Working Note 165 for fur‐
315 ther details and extra cautions. NPARAMS (input) INTE‐
316 GER Specifies the number of parameters set in PARAMS.
317 If .LE. 0, the PARAMS array is never referenced and
318 default values are used.
319
320 PARAMS (input / output) REAL array, dimension NPARAMS
321 Specifies algorithm parameters. If an entry is .LT. 0.0, then
322 that entry will be filled with default value used for that
323 parameter. Only positions up to NPARAMS are accessed; defaults
324 are used for higher-numbered parameters.
325 PARAMS(LA_LINRX_ITREF_I = 1) : Whether to perform iterative
326 refinement or not. Default: 1.0
327 = 0.0 : No refinement is performed, and no error bounds are
328 computed. = 1.0 : Use the double-precision refinement algo‐
329 rithm, possibly with doubled-single computations if the compi‐
330 lation environment does not support DOUBLE PRECISION. (other
331 values are reserved for future use) PARAMS(LA_LINRX_ITHRESH_I =
332 2) : Maximum number of residual computations allowed for
333 refinement. Default: 10
334 Aggressive: Set to 100 to permit convergence using approximate
335 factorizations or factorizations other than LU. If the factor‐
336 ization uses a technique other than Gaussian elimination, the
337 guarantees in err_bnds_norm and err_bnds_comp may no longer be
338 trustworthy. PARAMS(LA_LINRX_CWISE_I = 3) : Flag determining
339 if the code will attempt to find a solution with small compo‐
340 nentwise relative error in the double-precision algorithm.
341 Positive is true, 0.0 is false. Default: 1.0 (attempt compo‐
342 nentwise convergence)
343
344 WORK (workspace) REAL array, dimension (4*N)
345
346 IWORK (workspace) INTEGER array, dimension (N)
347
348 INFO (output) INTEGER
349 = 0: Successful exit. The solution to every right-hand side is
350 guaranteed. < 0: If INFO = -i, the i-th argument had an ille‐
351 gal value
352 > 0 and <= N: U(INFO,INFO) is exactly zero. The factorization
353 has been completed, but the factor U is exactly singular, so
354 the solution and error bounds could not be computed. RCOND = 0
355 is returned. = N+J: The solution corresponding to the Jth
356 right-hand side is not guaranteed. The solutions corresponding
357 to other right- hand sides K with K > J may not be guaranteed
358 as well, but only the first such right-hand side is reported.
359 If a small componentwise error is not requested (PARAMS(3) =
360 0.0) then the Jth right-hand side is the first with a normwise
361 error bound that is not guaranteed (the smallest J such that
362 ERR_BNDS_NORM(J,1) = 0.0). By default (PARAMS(3) = 1.0) the Jth
363 right-hand side is the first with either a normwise or compo‐
364 nentwise error bound that is not guaranteed (the smallest J
365 such that either ERR_BNDS_NORM(J,1) = 0.0 or ERR_BNDS_COMP(J,1)
366 = 0.0). See the definition of ERR_BNDS_NORM(:,1) and
367 ERR_BNDS_COMP(:,1). To get information about all of the right-
368 hand sides check ERR_BNDS_NORM or ERR_BNDS_COMP.
369
370
371
372 LAPACK driver routine (versionNo3v.e2m)ber 2008 SGBSVXX(1)