LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
dgbsv.f
Go to the documentation of this file.
1*> \brief <b> DGBSV computes the solution to system of linear equations A * X = B for GB matrices</b> (simple driver)
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download DGBSV + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dgbsv.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dgbsv.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dgbsv.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE DGBSV( N, KL, KU, NRHS, AB, LDAB, IPIV, B, LDB, INFO )
22*
23* .. Scalar Arguments ..
24* INTEGER INFO, KL, KU, LDAB, LDB, N, NRHS
25* ..
26* .. Array Arguments ..
27* INTEGER IPIV( * )
28* DOUBLE PRECISION AB( LDAB, * ), B( LDB, * )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> DGBSV computes the solution to a real system of linear equations
38*> A * X = B, where A is a band matrix of order N with KL subdiagonals
39*> and KU superdiagonals, and X and B are N-by-NRHS matrices.
40*>
41*> The LU decomposition with partial pivoting and row interchanges is
42*> used to factor A as A = L * U, where L is a product of permutation
43*> and unit lower triangular matrices with KL subdiagonals, and U is
44*> upper triangular with KL+KU superdiagonals. The factored form of A
45*> is then used to solve the system of equations A * X = B.
46*> \endverbatim
47*
48* Arguments:
49* ==========
50*
51*> \param[in] N
52*> \verbatim
53*> N is INTEGER
54*> The number of linear equations, i.e., the order of the
55*> matrix A. N >= 0.
56*> \endverbatim
57*>
58*> \param[in] KL
59*> \verbatim
60*> KL is INTEGER
61*> The number of subdiagonals within the band of A. KL >= 0.
62*> \endverbatim
63*>
64*> \param[in] KU
65*> \verbatim
66*> KU is INTEGER
67*> The number of superdiagonals within the band of A. KU >= 0.
68*> \endverbatim
69*>
70*> \param[in] NRHS
71*> \verbatim
72*> NRHS is INTEGER
73*> The number of right hand sides, i.e., the number of columns
74*> of the matrix B. NRHS >= 0.
75*> \endverbatim
76*>
77*> \param[in,out] AB
78*> \verbatim
79*> AB is DOUBLE PRECISION array, dimension (LDAB,N)
80*> On entry, the matrix A in band storage, in rows KL+1 to
81*> 2*KL+KU+1; rows 1 to KL of the array need not be set.
82*> The j-th column of A is stored in the j-th column of the
83*> array AB as follows:
84*> AB(KL+KU+1+i-j,j) = A(i,j) for max(1,j-KU)<=i<=min(N,j+KL)
85*> On exit, details of the factorization: U is stored as an
86*> upper triangular band matrix with KL+KU superdiagonals in
87*> rows 1 to KL+KU+1, and the multipliers used during the
88*> factorization are stored in rows KL+KU+2 to 2*KL+KU+1.
89*> See below for further details.
90*> \endverbatim
91*>
92*> \param[in] LDAB
93*> \verbatim
94*> LDAB is INTEGER
95*> The leading dimension of the array AB. LDAB >= 2*KL+KU+1.
96*> \endverbatim
97*>
98*> \param[out] IPIV
99*> \verbatim
100*> IPIV is INTEGER array, dimension (N)
101*> The pivot indices that define the permutation matrix P;
102*> row i of the matrix was interchanged with row IPIV(i).
103*> \endverbatim
104*>
105*> \param[in,out] B
106*> \verbatim
107*> B is DOUBLE PRECISION array, dimension (LDB,NRHS)
108*> On entry, the N-by-NRHS right hand side matrix B.
109*> On exit, if INFO = 0, the N-by-NRHS solution matrix X.
110*> \endverbatim
111*>
112*> \param[in] LDB
113*> \verbatim
114*> LDB is INTEGER
115*> The leading dimension of the array B. LDB >= max(1,N).
116*> \endverbatim
117*>
118*> \param[out] INFO
119*> \verbatim
120*> INFO is INTEGER
121*> = 0: successful exit
122*> < 0: if INFO = -i, the i-th argument had an illegal value
123*> > 0: if INFO = i, U(i,i) is exactly zero. The factorization
124*> has been completed, but the factor U is exactly
125*> singular, and the solution has not been computed.
126*> \endverbatim
127*
128* Authors:
129* ========
130*
131*> \author Univ. of Tennessee
132*> \author Univ. of California Berkeley
133*> \author Univ. of Colorado Denver
134*> \author NAG Ltd.
135*
136*> \ingroup gbsv
137*
138*> \par Further Details:
139* =====================
140*>
141*> \verbatim
142*>
143*> The band storage scheme is illustrated by the following example, when
144*> M = N = 6, KL = 2, KU = 1:
145*>
146*> On entry: On exit:
147*>
148*> * * * + + + * * * u14 u25 u36
149*> * * + + + + * * u13 u24 u35 u46
150*> * a12 a23 a34 a45 a56 * u12 u23 u34 u45 u56
151*> a11 a22 a33 a44 a55 a66 u11 u22 u33 u44 u55 u66
152*> a21 a32 a43 a54 a65 * m21 m32 m43 m54 m65 *
153*> a31 a42 a53 a64 * * m31 m42 m53 m64 * *
154*>
155*> Array elements marked * are not used by the routine; elements marked
156*> + need not be set on entry, but are required by the routine to store
157*> elements of U because of fill-in resulting from the row interchanges.
158*> \endverbatim
159*>
160* =====================================================================
161 SUBROUTINE dgbsv( N, KL, KU, NRHS, AB, LDAB, IPIV, B, LDB, INFO )
162*
163* -- LAPACK driver routine --
164* -- LAPACK is a software package provided by Univ. of Tennessee, --
165* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
166*
167* .. Scalar Arguments ..
168 INTEGER INFO, KL, KU, LDAB, LDB, N, NRHS
169* ..
170* .. Array Arguments ..
171 INTEGER IPIV( * )
172 DOUBLE PRECISION AB( LDAB, * ), B( LDB, * )
173* ..
174*
175* =====================================================================
176*
177* .. External Subroutines ..
178 EXTERNAL dgbtrf, dgbtrs, xerbla
179* ..
180* .. Intrinsic Functions ..
181 INTRINSIC max
182* ..
183* .. Executable Statements ..
184*
185* Test the input parameters.
186*
187 info = 0
188 IF( n.LT.0 ) THEN
189 info = -1
190 ELSE IF( kl.LT.0 ) THEN
191 info = -2
192 ELSE IF( ku.LT.0 ) THEN
193 info = -3
194 ELSE IF( nrhs.LT.0 ) THEN
195 info = -4
196 ELSE IF( ldab.LT.2*kl+ku+1 ) THEN
197 info = -6
198 ELSE IF( ldb.LT.max( n, 1 ) ) THEN
199 info = -9
200 END IF
201 IF( info.NE.0 ) THEN
202 CALL xerbla( 'DGBSV ', -info )
203 RETURN
204 END IF
205*
206* Compute the LU factorization of the band matrix A.
207*
208 CALL dgbtrf( n, n, kl, ku, ab, ldab, ipiv, info )
209 IF( info.EQ.0 ) THEN
210*
211* Solve the system A*X = B, overwriting B with X.
212*
213 CALL dgbtrs( 'No transpose', n, kl, ku, nrhs, ab, ldab, ipiv,
214 $ b, ldb, info )
215 END IF
216 RETURN
217*
218* End of DGBSV
219*
220 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgbsv(n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
DGBSV computes the solution to system of linear equations A * X = B for GB matrices (simple driver)
Definition dgbsv.f:162
subroutine dgbtrf(m, n, kl, ku, ab, ldab, ipiv, info)
DGBTRF
Definition dgbtrf.f:144
subroutine dgbtrs(trans, n, kl, ku, nrhs, ab, ldab, ipiv, b, ldb, info)
DGBTRS
Definition dgbtrs.f:138