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

◆ ztbmv()

subroutine ztbmv ( character  UPLO,
character  TRANS,
character  DIAG,
integer  N,
integer  K,
complex*16, dimension(lda,*)  A,
integer  LDA,
complex*16, dimension(*)  X,
integer  INCX 
)

ZTBMV

Purpose:
 ZTBMV  performs one of the matrix-vector operations

    x := A*x,   or   x := A**T*x,   or   x := A**H*x,

 where x is an n element vector and  A is an n by n unit, or non-unit,
 upper or lower triangular band matrix, with ( k + 1 ) diagonals.
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 operation to be performed as
           follows:

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

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

              TRANS = 'C' or 'c'   x := A**H*x.
[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 COMPLEX*16 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 COMPLEX*16 array, dimension at least
           ( 1 + ( n - 1 )*abs( INCX ) ).
           Before entry, the incremented array X must contain the n
           element vector x. On exit, X is overwritten with the
           transformed 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.
  The vector and matrix arguments are not referenced when N = 0, or M = 0

  -- 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 185 of file ztbmv.f.

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