LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
sgbtf2.f
Go to the documentation of this file.
1*> \brief \b SGBTF2 computes the LU factorization of a general band matrix using the unblocked version of the algorithm.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SGBTF2 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/sgbtf2.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/sgbtf2.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/sgbtf2.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SGBTF2( M, N, KL, KU, AB, LDAB, IPIV, INFO )
20*
21* .. Scalar Arguments ..
22* INTEGER INFO, KL, KU, LDAB, M, N
23* ..
24* .. Array Arguments ..
25* INTEGER IPIV( * )
26* REAL AB( LDAB, * )
27* ..
28*
29*
30*> \par Purpose:
31* =============
32*>
33*> \verbatim
34*>
35*> SGBTF2 computes an LU factorization of a real m-by-n band matrix A
36*> using partial pivoting with row interchanges.
37*>
38*> This is the unblocked version of the algorithm, calling Level 2 BLAS.
39*> \endverbatim
40*
41* Arguments:
42* ==========
43*
44*> \param[in] M
45*> \verbatim
46*> M is INTEGER
47*> The number of rows of the matrix A. M >= 0.
48*> \endverbatim
49*>
50*> \param[in] N
51*> \verbatim
52*> N is INTEGER
53*> The number of columns of the matrix A. N >= 0.
54*> \endverbatim
55*>
56*> \param[in] KL
57*> \verbatim
58*> KL is INTEGER
59*> The number of subdiagonals within the band of A. KL >= 0.
60*> \endverbatim
61*>
62*> \param[in] KU
63*> \verbatim
64*> KU is INTEGER
65*> The number of superdiagonals within the band of A. KU >= 0.
66*> \endverbatim
67*>
68*> \param[in,out] AB
69*> \verbatim
70*> AB is REAL array, dimension (LDAB,N)
71*> On entry, the matrix A in band storage, in rows KL+1 to
72*> 2*KL+KU+1; rows 1 to KL of the array need not be set.
73*> The j-th column of A is stored in the j-th column of the
74*> array AB as follows:
75*> AB(kl+ku+1+i-j,j) = A(i,j) for max(1,j-ku)<=i<=min(m,j+kl)
76*>
77*> On exit, details of the factorization: U is stored as an
78*> upper triangular band matrix with KL+KU superdiagonals in
79*> rows 1 to KL+KU+1, and the multipliers used during the
80*> factorization are stored in rows KL+KU+2 to 2*KL+KU+1.
81*> See below for further details.
82*> \endverbatim
83*>
84*> \param[in] LDAB
85*> \verbatim
86*> LDAB is INTEGER
87*> The leading dimension of the array AB. LDAB >= 2*KL+KU+1.
88*> \endverbatim
89*>
90*> \param[out] IPIV
91*> \verbatim
92*> IPIV is INTEGER array, dimension (min(M,N))
93*> The pivot indices; for 1 <= i <= min(M,N), row i of the
94*> matrix was interchanged with row IPIV(i).
95*> \endverbatim
96*>
97*> \param[out] INFO
98*> \verbatim
99*> INFO is INTEGER
100*> = 0: successful exit
101*> < 0: if INFO = -i, the i-th argument had an illegal value
102*> > 0: if INFO = +i, U(i,i) is exactly zero. The factorization
103*> has been completed, but the factor U is exactly
104*> singular, and division by zero will occur if it is used
105*> to solve a system of equations.
106*> \endverbatim
107*
108* Authors:
109* ========
110*
111*> \author Univ. of Tennessee
112*> \author Univ. of California Berkeley
113*> \author Univ. of Colorado Denver
114*> \author NAG Ltd.
115*
116*> \ingroup gbtf2
117*
118*> \par Further Details:
119* =====================
120*>
121*> \verbatim
122*>
123*> The band storage scheme is illustrated by the following example, when
124*> M = N = 6, KL = 2, KU = 1:
125*>
126*> On entry: On exit:
127*>
128*> * * * + + + * * * u14 u25 u36
129*> * * + + + + * * u13 u24 u35 u46
130*> * a12 a23 a34 a45 a56 * u12 u23 u34 u45 u56
131*> a11 a22 a33 a44 a55 a66 u11 u22 u33 u44 u55 u66
132*> a21 a32 a43 a54 a65 * m21 m32 m43 m54 m65 *
133*> a31 a42 a53 a64 * * m31 m42 m53 m64 * *
134*>
135*> Array elements marked * are not used by the routine; elements marked
136*> + need not be set on entry, but are required by the routine to store
137*> elements of U, because of fill-in resulting from the row
138*> interchanges.
139*> \endverbatim
140*>
141* =====================================================================
142 SUBROUTINE sgbtf2( M, N, KL, KU, AB, LDAB, IPIV, INFO )
143*
144* -- LAPACK computational routine --
145* -- LAPACK is a software package provided by Univ. of Tennessee, --
146* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
147*
148* .. Scalar Arguments ..
149 INTEGER INFO, KL, KU, LDAB, M, N
150* ..
151* .. Array Arguments ..
152 INTEGER IPIV( * )
153 REAL AB( LDAB, * )
154* ..
155*
156* =====================================================================
157*
158* .. Parameters ..
159 REAL ONE, ZERO
160 parameter( one = 1.0e+0, zero = 0.0e+0 )
161* ..
162* .. Local Scalars ..
163 INTEGER I, J, JP, JU, KM, KV
164* ..
165* .. External Functions ..
166 INTEGER ISAMAX
167 EXTERNAL isamax
168* ..
169* .. External Subroutines ..
170 EXTERNAL sger, sscal, sswap, xerbla
171* ..
172* .. Intrinsic Functions ..
173 INTRINSIC max, min
174* ..
175* .. Executable Statements ..
176*
177* KV is the number of superdiagonals in the factor U, allowing for
178* fill-in.
179*
180 kv = ku + kl
181*
182* Test the input parameters.
183*
184 info = 0
185 IF( m.LT.0 ) THEN
186 info = -1
187 ELSE IF( n.LT.0 ) THEN
188 info = -2
189 ELSE IF( kl.LT.0 ) THEN
190 info = -3
191 ELSE IF( ku.LT.0 ) THEN
192 info = -4
193 ELSE IF( ldab.LT.kl+kv+1 ) THEN
194 info = -6
195 END IF
196 IF( info.NE.0 ) THEN
197 CALL xerbla( 'SGBTF2', -info )
198 RETURN
199 END IF
200*
201* Quick return if possible
202*
203 IF( m.EQ.0 .OR. n.EQ.0 )
204 $ RETURN
205*
206* Gaussian elimination with partial pivoting
207*
208* Set fill-in elements in columns KU+2 to KV to zero.
209*
210 DO 20 j = ku + 2, min( kv, n )
211 DO 10 i = kv - j + 2, kl
212 ab( i, j ) = zero
213 10 CONTINUE
214 20 CONTINUE
215*
216* JU is the index of the last column affected by the current stage
217* of the factorization.
218*
219 ju = 1
220*
221 DO 40 j = 1, min( m, n )
222*
223* Set fill-in elements in column J+KV to zero.
224*
225 IF( j+kv.LE.n ) THEN
226 DO 30 i = 1, kl
227 ab( i, j+kv ) = zero
228 30 CONTINUE
229 END IF
230*
231* Find pivot and test for singularity. KM is the number of
232* subdiagonal elements in the current column.
233*
234 km = min( kl, m-j )
235 jp = isamax( km+1, ab( kv+1, j ), 1 )
236 ipiv( j ) = jp + j - 1
237 IF( ab( kv+jp, j ).NE.zero ) THEN
238 ju = max( ju, min( j+ku+jp-1, n ) )
239*
240* Apply interchange to columns J to JU.
241*
242 IF( jp.NE.1 )
243 $ CALL sswap( ju-j+1, ab( kv+jp, j ), ldab-1,
244 $ ab( kv+1, j ), ldab-1 )
245*
246 IF( km.GT.0 ) THEN
247*
248* Compute multipliers.
249*
250 CALL sscal( km, one / ab( kv+1, j ), ab( kv+2, j ),
251 $ 1 )
252*
253* Update trailing submatrix within the band.
254*
255 IF( ju.GT.j )
256 $ CALL sger( km, ju-j, -one, ab( kv+2, j ), 1,
257 $ ab( kv, j+1 ), ldab-1, ab( kv+1, j+1 ),
258 $ ldab-1 )
259 END IF
260 ELSE
261*
262* If pivot is zero, set INFO to the index of the pivot
263* unless a zero pivot has already been found.
264*
265 IF( info.EQ.0 )
266 $ info = j
267 END IF
268 40 CONTINUE
269 RETURN
270*
271* End of SGBTF2
272*
273 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sgbtf2(m, n, kl, ku, ab, ldab, ipiv, info)
SGBTF2 computes the LU factorization of a general band matrix using the unblocked version of the algo...
Definition sgbtf2.f:143
subroutine sger(m, n, alpha, x, incx, y, incy, a, lda)
SGER
Definition sger.f:130
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82