LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
cpbsvx.f
Go to the documentation of this file.
1*> \brief <b> CPBSVX computes the solution to system of linear equations A * X = B for OTHER matrices</b>
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download CPBSVX + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/cpbsvx.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/cpbsvx.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/cpbsvx.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE CPBSVX( FACT, UPLO, N, KD, NRHS, AB, LDAB, AFB, LDAFB,
20* EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR,
21* WORK, RWORK, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER EQUED, FACT, UPLO
25* INTEGER INFO, KD, LDAB, LDAFB, LDB, LDX, N, NRHS
26* REAL RCOND
27* ..
28* .. Array Arguments ..
29* REAL BERR( * ), FERR( * ), RWORK( * ), S( * )
30* COMPLEX AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ),
31* $ WORK( * ), X( LDX, * )
32* ..
33*
34*
35*> \par Purpose:
36* =============
37*>
38*> \verbatim
39*>
40*> CPBSVX uses the Cholesky factorization A = U**H*U or A = L*L**H to
41*> compute the solution to a complex system of linear equations
42*> A * X = B,
43*> where A is an N-by-N Hermitian positive definite band matrix and X
44*> and B are N-by-NRHS 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 = 'E', real scaling factors are computed to equilibrate
58*> the system:
59*> diag(S) * A * diag(S) * inv(diag(S)) * X = diag(S) * 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(S)*A*diag(S) and B by diag(S)*B.
63*>
64*> 2. If FACT = 'N' or 'E', the Cholesky decomposition is used to
65*> factor the matrix A (after equilibration if FACT = 'E') as
66*> A = U**H * U, if UPLO = 'U', or
67*> A = L * L**H, if UPLO = 'L',
68*> where U is an upper triangular band matrix, and L is a lower
69*> triangular band matrix.
70*>
71*> 3. If the leading principal minor of order i is not positive,
72*> then the routine returns with INFO = i. Otherwise, the factored
73*> form of A is used to estimate the condition number of the matrix
74*> A. If the reciprocal of the condition number is less than machine
75*> precision, INFO = N+1 is returned as a warning, but the routine
76*> still goes on to solve for X and compute error bounds as
77*> described below.
78*>
79*> 4. The system of equations is solved for X using the factored form
80*> of A.
81*>
82*> 5. Iterative refinement is applied to improve the computed solution
83*> matrix and calculate error bounds and backward error estimates
84*> for it.
85*>
86*> 6. If equilibration was used, the matrix X is premultiplied by
87*> diag(S) so that it solves the original system before
88*> equilibration.
89*> \endverbatim
90*
91* Arguments:
92* ==========
93*
94*> \param[in] FACT
95*> \verbatim
96*> FACT is CHARACTER*1
97*> Specifies whether or not the factored form of the matrix A is
98*> supplied on entry, and if not, whether the matrix A should be
99*> equilibrated before it is factored.
100*> = 'F': On entry, AFB contains the factored form of A.
101*> If EQUED = 'Y', the matrix A has been equilibrated
102*> with scaling factors given by S. AB and AFB will not
103*> be modified.
104*> = 'N': The matrix A will be copied to AFB and factored.
105*> = 'E': The matrix A will be equilibrated if necessary, then
106*> copied to AFB and factored.
107*> \endverbatim
108*>
109*> \param[in] UPLO
110*> \verbatim
111*> UPLO is CHARACTER*1
112*> = 'U': Upper triangle of A is stored;
113*> = 'L': Lower triangle of A is stored.
114*> \endverbatim
115*>
116*> \param[in] N
117*> \verbatim
118*> N is INTEGER
119*> The number of linear equations, i.e., the order of the
120*> matrix A. N >= 0.
121*> \endverbatim
122*>
123*> \param[in] KD
124*> \verbatim
125*> KD is INTEGER
126*> The number of superdiagonals of the matrix A if UPLO = 'U',
127*> or the number of subdiagonals if UPLO = 'L'. KD >= 0.
128*> \endverbatim
129*>
130*> \param[in] NRHS
131*> \verbatim
132*> NRHS is INTEGER
133*> The number of right-hand sides, i.e., the number of columns
134*> of the matrices B and X. NRHS >= 0.
135*> \endverbatim
136*>
137*> \param[in,out] AB
138*> \verbatim
139*> AB is COMPLEX array, dimension (LDAB,N)
140*> On entry, the upper or lower triangle of the Hermitian band
141*> matrix A, stored in the first KD+1 rows of the array, except
142*> if FACT = 'F' and EQUED = 'Y', then A must contain the
143*> equilibrated matrix diag(S)*A*diag(S). The j-th column of A
144*> is stored in the j-th column of the array AB as follows:
145*> if UPLO = 'U', AB(KD+1+i-j,j) = A(i,j) for max(1,j-KD)<=i<=j;
146*> if UPLO = 'L', AB(1+i-j,j) = A(i,j) for j<=i<=min(N,j+KD).
147*> See below for further details.
148*>
149*> On exit, if FACT = 'E' and EQUED = 'Y', A is overwritten by
150*> diag(S)*A*diag(S).
151*> \endverbatim
152*>
153*> \param[in] LDAB
154*> \verbatim
155*> LDAB is INTEGER
156*> The leading dimension of the array A. LDAB >= KD+1.
157*> \endverbatim
158*>
159*> \param[in,out] AFB
160*> \verbatim
161*> AFB is COMPLEX array, dimension (LDAFB,N)
162*> If FACT = 'F', then AFB is an input argument and on entry
163*> contains the triangular factor U or L from the Cholesky
164*> factorization A = U**H*U or A = L*L**H of the band matrix
165*> A, in the same storage format as A (see AB). If EQUED = 'Y',
166*> then AFB is the factored form of the equilibrated matrix A.
167*>
168*> If FACT = 'N', then AFB is an output argument and on exit
169*> returns the triangular factor U or L from the Cholesky
170*> factorization A = U**H*U or A = L*L**H.
171*>
172*> If FACT = 'E', then AFB is an output argument and on exit
173*> returns the triangular factor U or L from the Cholesky
174*> factorization A = U**H*U or A = L*L**H of the equilibrated
175*> matrix A (see the description of A for the form of the
176*> equilibrated matrix).
177*> \endverbatim
178*>
179*> \param[in] LDAFB
180*> \verbatim
181*> LDAFB is INTEGER
182*> The leading dimension of the array AFB. LDAFB >= KD+1.
183*> \endverbatim
184*>
185*> \param[in,out] EQUED
186*> \verbatim
187*> EQUED is CHARACTER*1
188*> Specifies the form of equilibration that was done.
189*> = 'N': No equilibration (always true if FACT = 'N').
190*> = 'Y': Equilibration was done, i.e., A has been replaced by
191*> diag(S) * A * diag(S).
192*> EQUED is an input argument if FACT = 'F'; otherwise, it is an
193*> output argument.
194*> \endverbatim
195*>
196*> \param[in,out] S
197*> \verbatim
198*> S is REAL array, dimension (N)
199*> The scale factors for A; not accessed if EQUED = 'N'. S is
200*> an input argument if FACT = 'F'; otherwise, S is an output
201*> argument. If FACT = 'F' and EQUED = 'Y', each element of S
202*> must be positive.
203*> \endverbatim
204*>
205*> \param[in,out] B
206*> \verbatim
207*> B is COMPLEX array, dimension (LDB,NRHS)
208*> On entry, the N-by-NRHS right hand side matrix B.
209*> On exit, if EQUED = 'N', B is not modified; if EQUED = 'Y',
210*> B is overwritten by diag(S) * B.
211*> \endverbatim
212*>
213*> \param[in] LDB
214*> \verbatim
215*> LDB is INTEGER
216*> The leading dimension of the array B. LDB >= max(1,N).
217*> \endverbatim
218*>
219*> \param[out] X
220*> \verbatim
221*> X is COMPLEX array, dimension (LDX,NRHS)
222*> If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X to
223*> the original system of equations. Note that if EQUED = 'Y',
224*> A and B are modified on exit, and the solution to the
225*> equilibrated system is inv(diag(S))*X.
226*> \endverbatim
227*>
228*> \param[in] LDX
229*> \verbatim
230*> LDX is INTEGER
231*> The leading dimension of the array X. LDX >= max(1,N).
232*> \endverbatim
233*>
234*> \param[out] RCOND
235*> \verbatim
236*> RCOND is REAL
237*> The estimate of the reciprocal condition number of the matrix
238*> A after equilibration (if done). If RCOND is less than the
239*> machine precision (in particular, if RCOND = 0), the matrix
240*> is singular to working precision. This condition is
241*> indicated by a return code of INFO > 0.
242*> \endverbatim
243*>
244*> \param[out] FERR
245*> \verbatim
246*> FERR is REAL array, dimension (NRHS)
247*> The estimated forward error bound for each solution vector
248*> X(j) (the j-th column of the solution matrix X).
249*> If XTRUE is the true solution corresponding to X(j), FERR(j)
250*> is an estimated upper bound for the magnitude of the largest
251*> element in (X(j) - XTRUE) divided by the magnitude of the
252*> largest element in X(j). The estimate is as reliable as
253*> the estimate for RCOND, and is almost always a slight
254*> overestimate of the true error.
255*> \endverbatim
256*>
257*> \param[out] BERR
258*> \verbatim
259*> BERR is REAL array, dimension (NRHS)
260*> The componentwise relative backward error of each solution
261*> vector X(j) (i.e., the smallest relative change in
262*> any element of A or B that makes X(j) an exact solution).
263*> \endverbatim
264*>
265*> \param[out] WORK
266*> \verbatim
267*> WORK is COMPLEX array, dimension (2*N)
268*> \endverbatim
269*>
270*> \param[out] RWORK
271*> \verbatim
272*> RWORK is REAL array, dimension (N)
273*> \endverbatim
274*>
275*> \param[out] INFO
276*> \verbatim
277*> INFO is INTEGER
278*> = 0: successful exit
279*> < 0: if INFO = -i, the i-th argument had an illegal value
280*> > 0: if INFO = i, and i is
281*> <= N: the leading principal minor of order i of A
282*> is not positive, so the factorization could not
283*> be completed, and the solution has not been
284*> computed. RCOND = 0 is returned.
285*> = N+1: U is nonsingular, but RCOND is less than machine
286*> precision, meaning that the matrix is singular
287*> to working precision. Nevertheless, the
288*> solution and error bounds are computed because
289*> there are a number of situations where the
290*> computed solution can be more accurate than the
291*> value of RCOND would suggest.
292*> \endverbatim
293*
294* Authors:
295* ========
296*
297*> \author Univ. of Tennessee
298*> \author Univ. of California Berkeley
299*> \author Univ. of Colorado Denver
300*> \author NAG Ltd.
301*
302*> \ingroup pbsvx
303*
304*> \par Further Details:
305* =====================
306*>
307*> \verbatim
308*>
309*> The band storage scheme is illustrated by the following example, when
310*> N = 6, KD = 2, and UPLO = 'U':
311*>
312*> Two-dimensional storage of the Hermitian matrix A:
313*>
314*> a11 a12 a13
315*> a22 a23 a24
316*> a33 a34 a35
317*> a44 a45 a46
318*> a55 a56
319*> (aij=conjg(aji)) a66
320*>
321*> Band storage of the upper triangle of A:
322*>
323*> * * a13 a24 a35 a46
324*> * a12 a23 a34 a45 a56
325*> a11 a22 a33 a44 a55 a66
326*>
327*> Similarly, if UPLO = 'L' the format of A is as follows:
328*>
329*> a11 a22 a33 a44 a55 a66
330*> a21 a32 a43 a54 a65 *
331*> a31 a42 a53 a64 * *
332*>
333*> Array elements marked * are not used by the routine.
334*> \endverbatim
335*>
336* =====================================================================
337 SUBROUTINE cpbsvx( FACT, UPLO, N, KD, NRHS, AB, LDAB, AFB,
338 $ LDAFB,
339 $ EQUED, S, B, LDB, X, LDX, RCOND, FERR, BERR,
340 $ WORK, RWORK, INFO )
341*
342* -- LAPACK driver routine --
343* -- LAPACK is a software package provided by Univ. of Tennessee, --
344* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
345*
346* .. Scalar Arguments ..
347 CHARACTER EQUED, FACT, UPLO
348 INTEGER INFO, KD, LDAB, LDAFB, LDB, LDX, N, NRHS
349 REAL RCOND
350* ..
351* .. Array Arguments ..
352 REAL BERR( * ), FERR( * ), RWORK( * ), S( * )
353 COMPLEX AB( LDAB, * ), AFB( LDAFB, * ), B( LDB, * ),
354 $ WORK( * ), X( LDX, * )
355* ..
356*
357* =====================================================================
358*
359* .. Parameters ..
360 REAL ZERO, ONE
361 PARAMETER ( ZERO = 0.0e+0, one = 1.0e+0 )
362* ..
363* .. Local Scalars ..
364 LOGICAL EQUIL, NOFACT, RCEQU, UPPER
365 INTEGER I, INFEQU, J, J1, J2
366 REAL AMAX, ANORM, BIGNUM, SCOND, SMAX, SMIN, SMLNUM
367* ..
368* .. External Functions ..
369 LOGICAL LSAME
370 REAL CLANHB, SLAMCH
371 EXTERNAL LSAME, CLANHB, SLAMCH
372* ..
373* .. External Subroutines ..
374 EXTERNAL ccopy, clacpy, claqhb, cpbcon, cpbequ,
375 $ cpbrfs,
377* ..
378* .. Intrinsic Functions ..
379 INTRINSIC max, min
380* ..
381* .. Executable Statements ..
382*
383 info = 0
384 nofact = lsame( fact, 'N' )
385 equil = lsame( fact, 'E' )
386 upper = lsame( uplo, 'U' )
387 IF( nofact .OR. equil ) THEN
388 equed = 'N'
389 rcequ = .false.
390 ELSE
391 rcequ = lsame( equed, 'Y' )
392 smlnum = slamch( 'Safe minimum' )
393 bignum = one / smlnum
394 END IF
395*
396* Test the input parameters.
397*
398 IF( .NOT.nofact .AND.
399 $ .NOT.equil .AND.
400 $ .NOT.lsame( fact, 'F' ) )
401 $ THEN
402 info = -1
403 ELSE IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
404 info = -2
405 ELSE IF( n.LT.0 ) THEN
406 info = -3
407 ELSE IF( kd.LT.0 ) THEN
408 info = -4
409 ELSE IF( nrhs.LT.0 ) THEN
410 info = -5
411 ELSE IF( ldab.LT.kd+1 ) THEN
412 info = -7
413 ELSE IF( ldafb.LT.kd+1 ) THEN
414 info = -9
415 ELSE IF( lsame( fact, 'F' ) .AND. .NOT.
416 $ ( rcequ .OR. lsame( equed, 'N' ) ) ) THEN
417 info = -10
418 ELSE
419 IF( rcequ ) THEN
420 smin = bignum
421 smax = zero
422 DO 10 j = 1, n
423 smin = min( smin, s( j ) )
424 smax = max( smax, s( j ) )
425 10 CONTINUE
426 IF( smin.LE.zero ) THEN
427 info = -11
428 ELSE IF( n.GT.0 ) THEN
429 scond = max( smin, smlnum ) / min( smax, bignum )
430 ELSE
431 scond = one
432 END IF
433 END IF
434 IF( info.EQ.0 ) THEN
435 IF( ldb.LT.max( 1, n ) ) THEN
436 info = -13
437 ELSE IF( ldx.LT.max( 1, n ) ) THEN
438 info = -15
439 END IF
440 END IF
441 END IF
442*
443 IF( info.NE.0 ) THEN
444 CALL xerbla( 'CPBSVX', -info )
445 RETURN
446 END IF
447*
448 IF( equil ) THEN
449*
450* Compute row and column scalings to equilibrate the matrix A.
451*
452 CALL cpbequ( uplo, n, kd, ab, ldab, s, scond, amax, infequ )
453 IF( infequ.EQ.0 ) THEN
454*
455* Equilibrate the matrix.
456*
457 CALL claqhb( uplo, n, kd, ab, ldab, s, scond, amax,
458 $ equed )
459 rcequ = lsame( equed, 'Y' )
460 END IF
461 END IF
462*
463* Scale the right-hand side.
464*
465 IF( rcequ ) THEN
466 DO 30 j = 1, nrhs
467 DO 20 i = 1, n
468 b( i, j ) = s( i )*b( i, j )
469 20 CONTINUE
470 30 CONTINUE
471 END IF
472*
473 IF( nofact .OR. equil ) THEN
474*
475* Compute the Cholesky factorization A = U**H *U or A = L*L**H.
476*
477 IF( upper ) THEN
478 DO 40 j = 1, n
479 j1 = max( j-kd, 1 )
480 CALL ccopy( j-j1+1, ab( kd+1-j+j1, j ), 1,
481 $ afb( kd+1-j+j1, j ), 1 )
482 40 CONTINUE
483 ELSE
484 DO 50 j = 1, n
485 j2 = min( j+kd, n )
486 CALL ccopy( j2-j+1, ab( 1, j ), 1, afb( 1, j ), 1 )
487 50 CONTINUE
488 END IF
489*
490 CALL cpbtrf( uplo, n, kd, afb, ldafb, info )
491*
492* Return if INFO is non-zero.
493*
494 IF( info.GT.0 )THEN
495 rcond = zero
496 RETURN
497 END IF
498 END IF
499*
500* Compute the norm of the matrix A.
501*
502 anorm = clanhb( '1', uplo, n, kd, ab, ldab, rwork )
503*
504* Compute the reciprocal of the condition number of A.
505*
506 CALL cpbcon( uplo, n, kd, afb, ldafb, anorm, rcond, work,
507 $ rwork,
508 $ info )
509*
510* Compute the solution matrix X.
511*
512 CALL clacpy( 'Full', n, nrhs, b, ldb, x, ldx )
513 CALL cpbtrs( uplo, n, kd, nrhs, afb, ldafb, x, ldx, info )
514*
515* Use iterative refinement to improve the computed solution and
516* compute error bounds and backward error estimates for it.
517*
518 CALL cpbrfs( uplo, n, kd, nrhs, ab, ldab, afb, ldafb, b, ldb,
519 $ x,
520 $ ldx, ferr, berr, work, rwork, info )
521*
522* Transform the solution matrix X to a solution of the original
523* system.
524*
525 IF( rcequ ) THEN
526 DO 70 j = 1, nrhs
527 DO 60 i = 1, n
528 x( i, j ) = s( i )*x( i, j )
529 60 CONTINUE
530 70 CONTINUE
531 DO 80 j = 1, nrhs
532 ferr( j ) = ferr( j ) / scond
533 80 CONTINUE
534 END IF
535*
536* Set INFO = N+1 if the matrix is singular to working precision.
537*
538 IF( rcond.LT.slamch( 'Epsilon' ) )
539 $ info = n + 1
540*
541 RETURN
542*
543* End of CPBSVX
544*
545 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ccopy(n, cx, incx, cy, incy)
CCOPY
Definition ccopy.f:81
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:101
subroutine claqhb(uplo, n, kd, ab, ldab, s, scond, amax, equed)
CLAQHB scales a Hermitian band matrix, using scaling factors computed by cpbequ.
Definition claqhb.f:140
subroutine cpbcon(uplo, n, kd, ab, ldab, anorm, rcond, work, rwork, info)
CPBCON
Definition cpbcon.f:131
subroutine cpbequ(uplo, n, kd, ab, ldab, s, scond, amax, info)
CPBEQU
Definition cpbequ.f:129
subroutine cpbrfs(uplo, n, kd, nrhs, ab, ldab, afb, ldafb, b, ldb, x, ldx, ferr, berr, work, rwork, info)
CPBRFS
Definition cpbrfs.f:187
subroutine cpbsvx(fact, uplo, n, kd, nrhs, ab, ldab, afb, ldafb, equed, s, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
CPBSVX computes the solution to system of linear equations A * X = B for OTHER matrices
Definition cpbsvx.f:341
subroutine cpbtrf(uplo, n, kd, ab, ldab, info)
CPBTRF
Definition cpbtrf.f:140
subroutine cpbtrs(uplo, n, kd, nrhs, ab, ldab, b, ldb, info)
CPBTRS
Definition cpbtrs.f:119