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