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

◆ ctbmv()

subroutine ctbmv ( character  uplo,
character  trans,
character  diag,
integer  n,
integer  k,
complex, dimension(lda,*)  a,
integer  lda,
complex, dimension(*)  x,
integer  incx 
)

CTBMV

Purpose:
 CTBMV  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 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 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 ctbmv.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 A(LDA,*),X(*)
197* ..
198*
199* =====================================================================
200*
201* .. Parameters ..
202 COMPLEX ZERO
203 parameter(zero= (0.0e+0,0.0e+0))
204* ..
205* .. Local Scalars ..
206 COMPLEX 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 conjg,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.
227 + .NOT.lsame(trans,'T') .AND.
228 + .NOT.lsame(trans,'C')) THEN
229 info = 2
230 ELSE IF (.NOT.lsame(diag,'U') .AND.
231 + .NOT.lsame(diag,'N')) THEN
232 info = 3
233 ELSE IF (n.LT.0) THEN
234 info = 4
235 ELSE IF (k.LT.0) THEN
236 info = 5
237 ELSE IF (lda.LT. (k+1)) THEN
238 info = 7
239 ELSE IF (incx.EQ.0) THEN
240 info = 9
241 END IF
242 IF (info.NE.0) THEN
243 CALL xerbla('CTBMV ',info)
244 RETURN
245 END IF
246*
247* Quick return if possible.
248*
249 IF (n.EQ.0) RETURN
250*
251 noconj = lsame(trans,'T')
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 sequentially with one pass through A.
265*
266 IF (lsame(trans,'N')) THEN
267*
268* Form x := A*x.
269*
270 IF (lsame(uplo,'U')) THEN
271 kplus1 = k + 1
272 IF (incx.EQ.1) THEN
273 DO 20 j = 1,n
274 IF (x(j).NE.zero) THEN
275 temp = x(j)
276 l = kplus1 - j
277 DO 10 i = max(1,j-k),j - 1
278 x(i) = x(i) + temp*a(l+i,j)
279 10 CONTINUE
280 IF (nounit) x(j) = x(j)*a(kplus1,j)
281 END IF
282 20 CONTINUE
283 ELSE
284 jx = kx
285 DO 40 j = 1,n
286 IF (x(jx).NE.zero) THEN
287 temp = x(jx)
288 ix = kx
289 l = kplus1 - j
290 DO 30 i = max(1,j-k),j - 1
291 x(ix) = x(ix) + temp*a(l+i,j)
292 ix = ix + incx
293 30 CONTINUE
294 IF (nounit) x(jx) = x(jx)*a(kplus1,j)
295 END IF
296 jx = jx + incx
297 IF (j.GT.k) kx = kx + incx
298 40 CONTINUE
299 END IF
300 ELSE
301 IF (incx.EQ.1) THEN
302 DO 60 j = n,1,-1
303 IF (x(j).NE.zero) THEN
304 temp = x(j)
305 l = 1 - j
306 DO 50 i = min(n,j+k),j + 1,-1
307 x(i) = x(i) + temp*a(l+i,j)
308 50 CONTINUE
309 IF (nounit) x(j) = x(j)*a(1,j)
310 END IF
311 60 CONTINUE
312 ELSE
313 kx = kx + (n-1)*incx
314 jx = kx
315 DO 80 j = n,1,-1
316 IF (x(jx).NE.zero) THEN
317 temp = x(jx)
318 ix = kx
319 l = 1 - j
320 DO 70 i = min(n,j+k),j + 1,-1
321 x(ix) = x(ix) + temp*a(l+i,j)
322 ix = ix - incx
323 70 CONTINUE
324 IF (nounit) x(jx) = x(jx)*a(1,j)
325 END IF
326 jx = jx - incx
327 IF ((n-j).GE.k) kx = kx - incx
328 80 CONTINUE
329 END IF
330 END IF
331 ELSE
332*
333* Form x := A**T*x or x := A**H*x.
334*
335 IF (lsame(uplo,'U')) THEN
336 kplus1 = k + 1
337 IF (incx.EQ.1) THEN
338 DO 110 j = n,1,-1
339 temp = x(j)
340 l = kplus1 - j
341 IF (noconj) THEN
342 IF (nounit) temp = temp*a(kplus1,j)
343 DO 90 i = j - 1,max(1,j-k),-1
344 temp = temp + a(l+i,j)*x(i)
345 90 CONTINUE
346 ELSE
347 IF (nounit) temp = temp*conjg(a(kplus1,j))
348 DO 100 i = j - 1,max(1,j-k),-1
349 temp = temp + conjg(a(l+i,j))*x(i)
350 100 CONTINUE
351 END IF
352 x(j) = temp
353 110 CONTINUE
354 ELSE
355 kx = kx + (n-1)*incx
356 jx = kx
357 DO 140 j = n,1,-1
358 temp = x(jx)
359 kx = kx - incx
360 ix = kx
361 l = kplus1 - j
362 IF (noconj) THEN
363 IF (nounit) temp = temp*a(kplus1,j)
364 DO 120 i = j - 1,max(1,j-k),-1
365 temp = temp + a(l+i,j)*x(ix)
366 ix = ix - incx
367 120 CONTINUE
368 ELSE
369 IF (nounit) temp = temp*conjg(a(kplus1,j))
370 DO 130 i = j - 1,max(1,j-k),-1
371 temp = temp + conjg(a(l+i,j))*x(ix)
372 ix = ix - incx
373 130 CONTINUE
374 END IF
375 x(jx) = temp
376 jx = jx - incx
377 140 CONTINUE
378 END IF
379 ELSE
380 IF (incx.EQ.1) THEN
381 DO 170 j = 1,n
382 temp = x(j)
383 l = 1 - j
384 IF (noconj) THEN
385 IF (nounit) temp = temp*a(1,j)
386 DO 150 i = j + 1,min(n,j+k)
387 temp = temp + a(l+i,j)*x(i)
388 150 CONTINUE
389 ELSE
390 IF (nounit) temp = temp*conjg(a(1,j))
391 DO 160 i = j + 1,min(n,j+k)
392 temp = temp + conjg(a(l+i,j))*x(i)
393 160 CONTINUE
394 END IF
395 x(j) = temp
396 170 CONTINUE
397 ELSE
398 jx = kx
399 DO 200 j = 1,n
400 temp = x(jx)
401 kx = kx + incx
402 ix = kx
403 l = 1 - j
404 IF (noconj) THEN
405 IF (nounit) temp = temp*a(1,j)
406 DO 180 i = j + 1,min(n,j+k)
407 temp = temp + a(l+i,j)*x(ix)
408 ix = ix + incx
409 180 CONTINUE
410 ELSE
411 IF (nounit) temp = temp*conjg(a(1,j))
412 DO 190 i = j + 1,min(n,j+k)
413 temp = temp + conjg(a(l+i,j))*x(ix)
414 ix = ix + incx
415 190 CONTINUE
416 END IF
417 x(jx) = temp
418 jx = jx + incx
419 200 CONTINUE
420 END IF
421 END IF
422 END IF
423*
424 RETURN
425*
426* End of CTBMV
427*
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: