LAPACK 3.12.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.
230 + .NOT.lsame(trans,'T') .AND.
231 + .NOT.lsame(trans,'C')) THEN
232 info = 2
233 ELSE IF (.NOT.lsame(diag,'U') .AND.
234 + .NOT.lsame(diag,'N')) THEN
235 info = 3
236 ELSE IF (n.LT.0) THEN
237 info = 4
238 ELSE IF (k.LT.0) THEN
239 info = 5
240 ELSE IF (lda.LT. (k+1)) THEN
241 info = 7
242 ELSE IF (incx.EQ.0) THEN
243 info = 9
244 END IF
245 IF (info.NE.0) THEN
246 CALL xerbla('DTBSV ',info)
247 RETURN
248 END IF
249*
250* Quick return if possible.
251*
252 IF (n.EQ.0) RETURN
253*
254 nounit = lsame(diag,'N')
255*
256* Set up the start point in X if the increment is not unity. This
257* will be ( N - 1 )*INCX too small for descending loops.
258*
259 IF (incx.LE.0) THEN
260 kx = 1 - (n-1)*incx
261 ELSE IF (incx.NE.1) THEN
262 kx = 1
263 END IF
264*
265* Start the operations. In this version the elements of A are
266* accessed by sequentially with one pass through A.
267*
268 IF (lsame(trans,'N')) THEN
269*
270* Form x := inv( A )*x.
271*
272 IF (lsame(uplo,'U')) THEN
273 kplus1 = k + 1
274 IF (incx.EQ.1) THEN
275 DO 20 j = n,1,-1
276 IF (x(j).NE.zero) THEN
277 l = kplus1 - j
278 IF (nounit) x(j) = x(j)/a(kplus1,j)
279 temp = x(j)
280 DO 10 i = j - 1,max(1,j-k),-1
281 x(i) = x(i) - temp*a(l+i,j)
282 10 CONTINUE
283 END IF
284 20 CONTINUE
285 ELSE
286 kx = kx + (n-1)*incx
287 jx = kx
288 DO 40 j = n,1,-1
289 kx = kx - incx
290 IF (x(jx).NE.zero) THEN
291 ix = kx
292 l = kplus1 - j
293 IF (nounit) x(jx) = x(jx)/a(kplus1,j)
294 temp = x(jx)
295 DO 30 i = j - 1,max(1,j-k),-1
296 x(ix) = x(ix) - temp*a(l+i,j)
297 ix = ix - incx
298 30 CONTINUE
299 END IF
300 jx = jx - incx
301 40 CONTINUE
302 END IF
303 ELSE
304 IF (incx.EQ.1) THEN
305 DO 60 j = 1,n
306 IF (x(j).NE.zero) THEN
307 l = 1 - j
308 IF (nounit) x(j) = x(j)/a(1,j)
309 temp = x(j)
310 DO 50 i = j + 1,min(n,j+k)
311 x(i) = x(i) - temp*a(l+i,j)
312 50 CONTINUE
313 END IF
314 60 CONTINUE
315 ELSE
316 jx = kx
317 DO 80 j = 1,n
318 kx = kx + incx
319 IF (x(jx).NE.zero) THEN
320 ix = kx
321 l = 1 - j
322 IF (nounit) x(jx) = x(jx)/a(1,j)
323 temp = x(jx)
324 DO 70 i = j + 1,min(n,j+k)
325 x(ix) = x(ix) - temp*a(l+i,j)
326 ix = ix + incx
327 70 CONTINUE
328 END IF
329 jx = jx + incx
330 80 CONTINUE
331 END IF
332 END IF
333 ELSE
334*
335* Form x := inv( A**T)*x.
336*
337 IF (lsame(uplo,'U')) THEN
338 kplus1 = k + 1
339 IF (incx.EQ.1) THEN
340 DO 100 j = 1,n
341 temp = x(j)
342 l = kplus1 - j
343 DO 90 i = max(1,j-k),j - 1
344 temp = temp - a(l+i,j)*x(i)
345 90 CONTINUE
346 IF (nounit) temp = temp/a(kplus1,j)
347 x(j) = temp
348 100 CONTINUE
349 ELSE
350 jx = kx
351 DO 120 j = 1,n
352 temp = x(jx)
353 ix = kx
354 l = kplus1 - j
355 DO 110 i = max(1,j-k),j - 1
356 temp = temp - a(l+i,j)*x(ix)
357 ix = ix + incx
358 110 CONTINUE
359 IF (nounit) temp = temp/a(kplus1,j)
360 x(jx) = temp
361 jx = jx + incx
362 IF (j.GT.k) kx = kx + incx
363 120 CONTINUE
364 END IF
365 ELSE
366 IF (incx.EQ.1) THEN
367 DO 140 j = n,1,-1
368 temp = x(j)
369 l = 1 - j
370 DO 130 i = min(n,j+k),j + 1,-1
371 temp = temp - a(l+i,j)*x(i)
372 130 CONTINUE
373 IF (nounit) temp = temp/a(1,j)
374 x(j) = temp
375 140 CONTINUE
376 ELSE
377 kx = kx + (n-1)*incx
378 jx = kx
379 DO 160 j = n,1,-1
380 temp = x(jx)
381 ix = kx
382 l = 1 - j
383 DO 150 i = min(n,j+k),j + 1,-1
384 temp = temp - a(l+i,j)*x(ix)
385 ix = ix - incx
386 150 CONTINUE
387 IF (nounit) temp = temp/a(1,j)
388 x(jx) = temp
389 jx = jx - incx
390 IF ((n-j).GE.k) kx = kx - incx
391 160 CONTINUE
392 END IF
393 END IF
394 END IF
395*
396 RETURN
397*
398* End of DTBSV
399*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: