LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine sgbtf2 ( integer  M,
integer  N,
integer  KL,
integer  KU,
real, dimension( ldab, * )  AB,
integer  LDAB,
integer, dimension( * )  IPIV,
integer  INFO 
)

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

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

Purpose:
 SGBTF2 computes an LU factorization of a real 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 REAL 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.
Date
September 2012
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 147 of file sgbtf2.f.

147 *
148 * -- LAPACK computational routine (version 3.4.2) --
149 * -- LAPACK is a software package provided by Univ. of Tennessee, --
150 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
151 * September 2012
152 *
153 * .. Scalar Arguments ..
154  INTEGER info, kl, ku, ldab, m, n
155 * ..
156 * .. Array Arguments ..
157  INTEGER ipiv( * )
158  REAL ab( ldab, * )
159 * ..
160 *
161 * =====================================================================
162 *
163 * .. Parameters ..
164  REAL one, zero
165  parameter ( one = 1.0e+0, zero = 0.0e+0 )
166 * ..
167 * .. Local Scalars ..
168  INTEGER i, j, jp, ju, km, kv
169 * ..
170 * .. External Functions ..
171  INTEGER isamax
172  EXTERNAL isamax
173 * ..
174 * .. External Subroutines ..
175  EXTERNAL sger, sscal, sswap, xerbla
176 * ..
177 * .. Intrinsic Functions ..
178  INTRINSIC max, min
179 * ..
180 * .. Executable Statements ..
181 *
182 * KV is the number of superdiagonals in the factor U, allowing for
183 * fill-in.
184 *
185  kv = ku + kl
186 *
187 * Test the input parameters.
188 *
189  info = 0
190  IF( m.LT.0 ) THEN
191  info = -1
192  ELSE IF( n.LT.0 ) THEN
193  info = -2
194  ELSE IF( kl.LT.0 ) THEN
195  info = -3
196  ELSE IF( ku.LT.0 ) THEN
197  info = -4
198  ELSE IF( ldab.LT.kl+kv+1 ) THEN
199  info = -6
200  END IF
201  IF( info.NE.0 ) THEN
202  CALL xerbla( 'SGBTF2', -info )
203  RETURN
204  END IF
205 *
206 * Quick return if possible
207 *
208  IF( m.EQ.0 .OR. n.EQ.0 )
209  $ RETURN
210 *
211 * Gaussian elimination with partial pivoting
212 *
213 * Set fill-in elements in columns KU+2 to KV to zero.
214 *
215  DO 20 j = ku + 2, min( kv, n )
216  DO 10 i = kv - j + 2, kl
217  ab( i, j ) = zero
218  10 CONTINUE
219  20 CONTINUE
220 *
221 * JU is the index of the last column affected by the current stage
222 * of the factorization.
223 *
224  ju = 1
225 *
226  DO 40 j = 1, min( m, n )
227 *
228 * Set fill-in elements in column J+KV to zero.
229 *
230  IF( j+kv.LE.n ) THEN
231  DO 30 i = 1, kl
232  ab( i, j+kv ) = zero
233  30 CONTINUE
234  END IF
235 *
236 * Find pivot and test for singularity. KM is the number of
237 * subdiagonal elements in the current column.
238 *
239  km = min( kl, m-j )
240  jp = isamax( km+1, ab( kv+1, j ), 1 )
241  ipiv( j ) = jp + j - 1
242  IF( ab( kv+jp, j ).NE.zero ) THEN
243  ju = max( ju, min( j+ku+jp-1, n ) )
244 *
245 * Apply interchange to columns J to JU.
246 *
247  IF( jp.NE.1 )
248  $ CALL sswap( ju-j+1, ab( kv+jp, j ), ldab-1,
249  $ ab( kv+1, j ), ldab-1 )
250 *
251  IF( km.GT.0 ) THEN
252 *
253 * Compute multipliers.
254 *
255  CALL sscal( km, one / ab( kv+1, j ), ab( kv+2, j ), 1 )
256 *
257 * Update trailing submatrix within the band.
258 *
259  IF( ju.GT.j )
260  $ CALL sger( km, ju-j, -one, ab( kv+2, j ), 1,
261  $ ab( kv, j+1 ), ldab-1, ab( kv+1, j+1 ),
262  $ ldab-1 )
263  END IF
264  ELSE
265 *
266 * If pivot is zero, set INFO to the index of the pivot
267 * unless a zero pivot has already been found.
268 *
269  IF( info.EQ.0 )
270  $ info = j
271  END IF
272  40 CONTINUE
273  RETURN
274 *
275 * End of SGBTF2
276 *
subroutine sger(M, N, ALPHA, X, INCX, Y, INCY, A, LDA)
SGER
Definition: sger.f:132
integer function isamax(N, SX, INCX)
ISAMAX
Definition: isamax.f:53
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
subroutine sswap(N, SX, INCX, SY, INCY)
SSWAP
Definition: sswap.f:53

Here is the call graph for this function:

Here is the caller graph for this function: