1CGBSVXX(1) LAPACK driver routine (version 3.2)                      CGBSVXX(1)
2
3
4

NAME

6       CGBSVXX - CGBSVXX use the LU factorization to compute the solution to a
7       complex system of linear equations A * X = B, where  A  is  an   N-by-N
8       matrix and X and B are N-by-NRHS matrices
9

SYNOPSIS

11       SUBROUTINE CGBSVXX( 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, RWORK, 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( * )
26
27           COMPLEX         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                           RWORK( * )
33

PURPOSE

35          CGBSVXX uses the LU factorization to compute the solution to a
36          complex system of linear equations  A * X = B,  where A is 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. CGBSVXX 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          CGBSVXX 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          CGBSVXX 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 CGBSVXX would itself produce.
51

DESCRIPTION

53          The following steps are performed:
54          1. If FACT = 'E', real scaling factors are computed to equilibrate
55          the system:
56            TRANS = 'N':  diag(R)*A*diag(C)     *inv(diag(C))*X = diag(R)*B
57            TRANS = 'T': (diag(R)*A*diag(C))**T *inv(diag(R))*X = diag(C)*B
58            TRANS = 'C': (diag(R)*A*diag(C))**H *inv(diag(R))*X = diag(C)*B
59          Whether or not the system will be equilibrated depends on the
60          scaling of the matrix A, but if equilibration is used, A is
61          overwritten by diag(R)*A*diag(C) and B by diag(R)*B (if TRANS='N')
62          or diag(C)*B (if TRANS = 'T' or 'C').
63          2. If FACT = 'N' or 'E', the LU decomposition is used to factor
64          the matrix A (after equilibration if FACT = 'E') as
65            A = P * L * U,
66          where P is a permutation matrix, L is a unit lower triangular
67          matrix, and U is upper triangular.
68          3. If some U(i,i)=0, so that U is exactly singular, then the
69          routine returns with INFO = i. Otherwise, the factored form of A
70          is used to estimate the condition number of the matrix A (see
71          argument RCOND). If the reciprocal of the condition number is less
72          than machine precision, the routine still goes on to solve for X
73          and compute error bounds as described below.
74          4. The system of equations is solved for X using the factored form
75          of A.
76          5. By default (unless PARAMS(LA_LINRX_ITREF_I) is set to zero),
77          the routine will use iterative refinement to try to get a small
78          error and error bounds.  Refinement calculates the residual to at
79          least twice the working precision.
80          6. If equilibration was used, the matrix X is premultiplied by
81          diag(C) (if TRANS = 'N') or diag(R) (if TRANS = 'T' or 'C') so
82          that it solves the original system before equilibration.
83

ARGUMENTS

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