LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ dtbsv()

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, 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, 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.
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 188 of file dtbsv.f.

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