LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages

◆ zgbtf2()

subroutine zgbtf2 ( integer m,
integer n,
integer kl,
integer ku,
complex*16, dimension( ldab, * ) ab,
integer ldab,
integer, dimension( * ) ipiv,
integer info )

ZGBTF2 computes the LU factorization of a general band matrix using the unblocked version of the algorithm.

Download ZGBTF2 + dependencies [TGZ] [ZIP] [TXT]

Purpose:
!> !> ZGBTF2 computes an LU factorization of a complex m-by-n band matrix !> A using partial pivoting with row interchanges. !> !> This is the unblocked version of the algorithm, calling Level 2 BLAS. !>
Parameters
[in]M
!> M is INTEGER !> The number of rows of the matrix A. M >= 0. !>
[in]N
!> N is INTEGER !> The number of columns of the matrix A. N >= 0. !>
[in]KL
!> KL is INTEGER !> The number of subdiagonals within the band of A. KL >= 0. !>
[in]KU
!> KU is INTEGER !> The number of superdiagonals within the band of A. KU >= 0. !>
[in,out]AB
!> AB is COMPLEX*16 array, dimension (LDAB,N) !> On entry, the matrix A in band storage, in rows KL+1 to !> 2*KL+KU+1; rows 1 to KL of the array need not be set. !> The j-th column of A is stored in the j-th column of the !> array AB as follows: !> AB(kl+ku+1+i-j,j) = A(i,j) for max(1,j-ku)<=i<=min(m,j+kl) !> !> On exit, details of the factorization: U is stored as an !> upper triangular band matrix with KL+KU superdiagonals in !> rows 1 to KL+KU+1, and the multipliers used during the !> factorization are stored in rows KL+KU+2 to 2*KL+KU+1. !> See below for further details. !>
[in]LDAB
!> LDAB is INTEGER !> The leading dimension of the array AB. LDAB >= 2*KL+KU+1. !>
[out]IPIV
!> IPIV is INTEGER array, dimension (min(M,N)) !> The pivot indices; for 1 <= i <= min(M,N), row i of the !> matrix was interchanged with row IPIV(i). !>
[out]INFO
!> INFO is INTEGER !> = 0: successful exit !> < 0: if INFO = -i, the i-th argument had an illegal value !> > 0: if INFO = +i, U(i,i) is exactly zero. The factorization !> has been completed, but the factor U is exactly !> singular, and division by zero will occur if it is used !> to solve a system of equations. !>
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!> !> The band storage scheme is illustrated by the following example, when !> M = N = 6, KL = 2, KU = 1: !> !> On entry: On exit: !> !> * * * + + + * * * u14 u25 u36 !> * * + + + + * * u13 u24 u35 u46 !> * a12 a23 a34 a45 a56 * u12 u23 u34 u45 u56 !> a11 a22 a33 a44 a55 a66 u11 u22 u33 u44 u55 u66 !> a21 a32 a43 a54 a65 * m21 m32 m43 m54 m65 * !> a31 a42 a53 a64 * * m31 m42 m53 m64 * * !> !> Array elements marked * are not used by the routine; elements marked !> + need not be set on entry, but are required by the routine to store !> elements of U, because of fill-in resulting from the row !> interchanges. !>

Definition at line 142 of file zgbtf2.f.

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 COMPLEX*16 AB( LDAB, * )
154* ..
155*
156* =====================================================================
157*
158* .. Parameters ..
159 COMPLEX*16 ONE, ZERO
160 parameter( one = ( 1.0d+0, 0.0d+0 ),
161 $ zero = ( 0.0d+0, 0.0d+0 ) )
162* ..
163* .. Local Scalars ..
164 INTEGER I, J, JP, JU, KM, KV
165* ..
166* .. External Functions ..
167 INTEGER IZAMAX
168 EXTERNAL izamax
169* ..
170* .. External Subroutines ..
171 EXTERNAL xerbla, zgeru, zscal, zswap
172* ..
173* .. Intrinsic Functions ..
174 INTRINSIC max, min
175* ..
176* .. Executable Statements ..
177*
178* KV is the number of superdiagonals in the factor U, allowing for
179* fill-in.
180*
181 kv = ku + kl
182*
183* Test the input parameters.
184*
185 info = 0
186 IF( m.LT.0 ) THEN
187 info = -1
188 ELSE IF( n.LT.0 ) THEN
189 info = -2
190 ELSE IF( kl.LT.0 ) THEN
191 info = -3
192 ELSE IF( ku.LT.0 ) THEN
193 info = -4
194 ELSE IF( ldab.LT.kl+kv+1 ) THEN
195 info = -6
196 END IF
197 IF( info.NE.0 ) THEN
198 CALL xerbla( 'ZGBTF2', -info )
199 RETURN
200 END IF
201*
202* Quick return if possible
203*
204 IF( m.EQ.0 .OR. n.EQ.0 )
205 $ RETURN
206*
207* Gaussian elimination with partial pivoting
208*
209* Set fill-in elements in columns KU+2 to KV to zero.
210*
211 DO 20 j = ku + 2, min( kv, n )
212 DO 10 i = kv - j + 2, kl
213 ab( i, j ) = zero
214 10 CONTINUE
215 20 CONTINUE
216*
217* JU is the index of the last column affected by the current stage
218* of the factorization.
219*
220 ju = 1
221*
222 DO 40 j = 1, min( m, n )
223*
224* Set fill-in elements in column J+KV to zero.
225*
226 IF( j+kv.LE.n ) THEN
227 DO 30 i = 1, kl
228 ab( i, j+kv ) = zero
229 30 CONTINUE
230 END IF
231*
232* Find pivot and test for singularity. KM is the number of
233* subdiagonal elements in the current column.
234*
235 km = min( kl, m-j )
236 jp = izamax( km+1, ab( kv+1, j ), 1 )
237 ipiv( j ) = jp + j - 1
238 IF( ab( kv+jp, j ).NE.zero ) THEN
239 ju = max( ju, min( j+ku+jp-1, n ) )
240*
241* Apply interchange to columns J to JU.
242*
243 IF( jp.NE.1 )
244 $ CALL zswap( ju-j+1, ab( kv+jp, j ), ldab-1,
245 $ ab( kv+1, j ), ldab-1 )
246 IF( km.GT.0 ) THEN
247*
248* Compute multipliers.
249*
250 CALL zscal( 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 zgeru( 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 ZGBTF2
272*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgeru(m, n, alpha, x, incx, y, incy, a, lda)
ZGERU
Definition zgeru.f:130
integer function izamax(n, zx, incx)
IZAMAX
Definition izamax.f:71
subroutine zscal(n, za, zx, incx)
ZSCAL
Definition zscal.f:78
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81
Here is the call graph for this function:
Here is the caller graph for this function: