LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dtbsv ( character  UPLO,
character  TRANS,
character  DIAG,
integer  N,
integer  K,
double precision, dimension(lda,*)  A,
integer  LDA,
double precision, dimension(*)  X,
integer  INCX 
)

DTBSV

Purpose:
 DTBSV  solves one of the systems of equations

    A*x = b,   or   A**T*x = b,

 where b and x are n element vectors and A is an n by n unit, or
 non-unit, upper or lower triangular band matrix, with ( k + 1 )
 diagonals.

 No test for singularity or near-singularity is included in this
 routine. Such tests must be performed before calling this routine.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
           On entry, UPLO specifies whether the matrix is an upper or
           lower triangular matrix as follows:

              UPLO = 'U' or 'u'   A is an upper triangular matrix.

              UPLO = 'L' or 'l'   A is a lower triangular matrix.
[in]TRANS
          TRANS is CHARACTER*1
           On entry, TRANS specifies the equations to be solved as
           follows:

              TRANS = 'N' or 'n'   A*x = b.

              TRANS = 'T' or 't'   A**T*x = b.

              TRANS = 'C' or 'c'   A**T*x = b.
[in]DIAG
          DIAG is CHARACTER*1
           On entry, DIAG specifies whether or not A is unit
           triangular as follows:

              DIAG = 'U' or 'u'   A is assumed to be unit triangular.

              DIAG = 'N' or 'n'   A is not assumed to be unit
                                  triangular.
[in]N
          N is INTEGER
           On entry, N specifies the order of the matrix A.
           N must be at least zero.
[in]K
          K is INTEGER
           On entry with UPLO = 'U' or 'u', K specifies the number of
           super-diagonals of the matrix A.
           On entry with UPLO = 'L' or 'l', K specifies the number of
           sub-diagonals of the matrix A.
           K must satisfy  0 .le. K.
[in]A
          A is DOUBLE PRECISION array of DIMENSION ( LDA, n ).
           Before entry with UPLO = 'U' or 'u', the leading ( k + 1 )
           by n part of the array A must contain the upper triangular
           band part of the matrix of coefficients, supplied column by
           column, with the leading diagonal of the matrix in row
           ( k + 1 ) of the array, the first super-diagonal starting at
           position 2 in row k, and so on. The top left k by k triangle
           of the array A is not referenced.
           The following program segment will transfer an upper
           triangular band matrix from conventional full matrix storage
           to band storage:

                 DO 20, J = 1, N
                    M = K + 1 - J
                    DO 10, I = MAX( 1, J - K ), J
                       A( M + I, J ) = matrix( I, J )
              10    CONTINUE
              20 CONTINUE

           Before entry with UPLO = 'L' or 'l', the leading ( k + 1 )
           by n part of the array A must contain the lower triangular
           band part of the matrix of coefficients, supplied column by
           column, with the leading diagonal of the matrix in row 1 of
           the array, the first sub-diagonal starting at position 1 in
           row 2, and so on. The bottom right k by k triangle of the
           array A is not referenced.
           The following program segment will transfer a lower
           triangular band matrix from conventional full matrix storage
           to band storage:

                 DO 20, J = 1, N
                    M = 1 - J
                    DO 10, I = J, MIN( N, J + K )
                       A( M + I, J ) = matrix( I, J )
              10    CONTINUE
              20 CONTINUE

           Note that when DIAG = 'U' or 'u' the elements of the array A
           corresponding to the diagonal elements of the matrix are not
           referenced, but are assumed to be unity.
[in]LDA
          LDA is INTEGER
           On entry, LDA specifies the first dimension of A as declared
           in the calling (sub) program. LDA must be at least
           ( k + 1 ).
[in,out]X
          X is DOUBLE PRECISION array of dimension at least
           ( 1 + ( n - 1 )*abs( INCX ) ).
           Before entry, the incremented array X must contain the n
           element right-hand side vector b. On exit, X is overwritten
           with the solution vector x.
[in]INCX
          INCX is INTEGER
           On entry, INCX specifies the increment for the elements of
           X. INCX must not be zero.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011
Further Details:
  Level 2 Blas routine.

  -- Written on 22-October-1986.
     Jack Dongarra, Argonne National Lab.
     Jeremy Du Croz, Nag Central Office.
     Sven Hammarling, Nag Central Office.
     Richard Hanson, Sandia National Labs.

Definition at line 191 of file dtbsv.f.

191 *
192 * -- Reference BLAS level2 routine (version 3.4.0) --
193 * -- Reference BLAS is a software package provided by Univ. of Tennessee, --
194 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
195 * November 2011
196 *
197 * .. Scalar Arguments ..
198  INTEGER incx,k,lda,n
199  CHARACTER diag,trans,uplo
200 * ..
201 * .. Array Arguments ..
202  DOUBLE PRECISION a(lda,*),x(*)
203 * ..
204 *
205 * =====================================================================
206 *
207 * .. Parameters ..
208  DOUBLE PRECISION zero
209  parameter(zero=0.0d+0)
210 * ..
211 * .. Local Scalars ..
212  DOUBLE PRECISION temp
213  INTEGER i,info,ix,j,jx,kplus1,kx,l
214  LOGICAL nounit
215 * ..
216 * .. External Functions ..
217  LOGICAL lsame
218  EXTERNAL lsame
219 * ..
220 * .. External Subroutines ..
221  EXTERNAL xerbla
222 * ..
223 * .. Intrinsic Functions ..
224  INTRINSIC max,min
225 * ..
226 *
227 * Test the input parameters.
228 *
229  info = 0
230  IF (.NOT.lsame(uplo,'U') .AND. .NOT.lsame(uplo,'L')) THEN
231  info = 1
232  ELSE IF (.NOT.lsame(trans,'N') .AND. .NOT.lsame(trans,'T') .AND.
233  + .NOT.lsame(trans,'C')) THEN
234  info = 2
235  ELSE IF (.NOT.lsame(diag,'U') .AND. .NOT.lsame(diag,'N')) THEN
236  info = 3
237  ELSE IF (n.LT.0) THEN
238  info = 4
239  ELSE IF (k.LT.0) THEN
240  info = 5
241  ELSE IF (lda.LT. (k+1)) THEN
242  info = 7
243  ELSE IF (incx.EQ.0) THEN
244  info = 9
245  END IF
246  IF (info.NE.0) THEN
247  CALL xerbla('DTBSV ',info)
248  RETURN
249  END IF
250 *
251 * Quick return if possible.
252 *
253  IF (n.EQ.0) RETURN
254 *
255  nounit = lsame(diag,'N')
256 *
257 * Set up the start point in X if the increment is not unity. This
258 * will be ( N - 1 )*INCX too small for descending loops.
259 *
260  IF (incx.LE.0) THEN
261  kx = 1 - (n-1)*incx
262  ELSE IF (incx.NE.1) THEN
263  kx = 1
264  END IF
265 *
266 * Start the operations. In this version the elements of A are
267 * accessed by sequentially with one pass through A.
268 *
269  IF (lsame(trans,'N')) THEN
270 *
271 * Form x := inv( A )*x.
272 *
273  IF (lsame(uplo,'U')) THEN
274  kplus1 = k + 1
275  IF (incx.EQ.1) THEN
276  DO 20 j = n,1,-1
277  IF (x(j).NE.zero) THEN
278  l = kplus1 - j
279  IF (nounit) x(j) = x(j)/a(kplus1,j)
280  temp = x(j)
281  DO 10 i = j - 1,max(1,j-k),-1
282  x(i) = x(i) - temp*a(l+i,j)
283  10 CONTINUE
284  END IF
285  20 CONTINUE
286  ELSE
287  kx = kx + (n-1)*incx
288  jx = kx
289  DO 40 j = n,1,-1
290  kx = kx - incx
291  IF (x(jx).NE.zero) THEN
292  ix = kx
293  l = kplus1 - j
294  IF (nounit) x(jx) = x(jx)/a(kplus1,j)
295  temp = x(jx)
296  DO 30 i = j - 1,max(1,j-k),-1
297  x(ix) = x(ix) - temp*a(l+i,j)
298  ix = ix - incx
299  30 CONTINUE
300  END IF
301  jx = jx - incx
302  40 CONTINUE
303  END IF
304  ELSE
305  IF (incx.EQ.1) THEN
306  DO 60 j = 1,n
307  IF (x(j).NE.zero) THEN
308  l = 1 - j
309  IF (nounit) x(j) = x(j)/a(1,j)
310  temp = x(j)
311  DO 50 i = j + 1,min(n,j+k)
312  x(i) = x(i) - temp*a(l+i,j)
313  50 CONTINUE
314  END IF
315  60 CONTINUE
316  ELSE
317  jx = kx
318  DO 80 j = 1,n
319  kx = kx + incx
320  IF (x(jx).NE.zero) THEN
321  ix = kx
322  l = 1 - j
323  IF (nounit) x(jx) = x(jx)/a(1,j)
324  temp = x(jx)
325  DO 70 i = j + 1,min(n,j+k)
326  x(ix) = x(ix) - temp*a(l+i,j)
327  ix = ix + incx
328  70 CONTINUE
329  END IF
330  jx = jx + incx
331  80 CONTINUE
332  END IF
333  END IF
334  ELSE
335 *
336 * Form x := inv( A**T)*x.
337 *
338  IF (lsame(uplo,'U')) THEN
339  kplus1 = k + 1
340  IF (incx.EQ.1) THEN
341  DO 100 j = 1,n
342  temp = x(j)
343  l = kplus1 - j
344  DO 90 i = max(1,j-k),j - 1
345  temp = temp - a(l+i,j)*x(i)
346  90 CONTINUE
347  IF (nounit) temp = temp/a(kplus1,j)
348  x(j) = temp
349  100 CONTINUE
350  ELSE
351  jx = kx
352  DO 120 j = 1,n
353  temp = x(jx)
354  ix = kx
355  l = kplus1 - j
356  DO 110 i = max(1,j-k),j - 1
357  temp = temp - a(l+i,j)*x(ix)
358  ix = ix + incx
359  110 CONTINUE
360  IF (nounit) temp = temp/a(kplus1,j)
361  x(jx) = temp
362  jx = jx + incx
363  IF (j.GT.k) kx = kx + incx
364  120 CONTINUE
365  END IF
366  ELSE
367  IF (incx.EQ.1) THEN
368  DO 140 j = n,1,-1
369  temp = x(j)
370  l = 1 - j
371  DO 130 i = min(n,j+k),j + 1,-1
372  temp = temp - a(l+i,j)*x(i)
373  130 CONTINUE
374  IF (nounit) temp = temp/a(1,j)
375  x(j) = temp
376  140 CONTINUE
377  ELSE
378  kx = kx + (n-1)*incx
379  jx = kx
380  DO 160 j = n,1,-1
381  temp = x(jx)
382  ix = kx
383  l = 1 - j
384  DO 150 i = min(n,j+k),j + 1,-1
385  temp = temp - a(l+i,j)*x(ix)
386  ix = ix - incx
387  150 CONTINUE
388  IF (nounit) temp = temp/a(1,j)
389  x(jx) = temp
390  jx = jx - incx
391  IF ((n-j).GE.k) kx = kx - incx
392  160 CONTINUE
393  END IF
394  END IF
395  END IF
396 *
397  RETURN
398 *
399 * End of DTBSV .
400 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: