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