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