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

◆ cgbmv()

subroutine cgbmv ( character trans,
integer m,
integer n,
integer kl,
integer ku,
complex alpha,
complex, dimension(lda,*) a,
integer lda,
complex, dimension(*) x,
integer incx,
complex beta,
complex, dimension(*) y,
integer incy )

CGBMV

Purpose:
!>
!> CGBMV  performs one of the matrix-vector operations
!>
!>    y := alpha*A*x + beta*y,   or   y := alpha*A**T*x + beta*y,   or
!>
!>    y := alpha*A**H*x + beta*y,
!>
!> where alpha and beta are scalars, x and y are vectors and A is an
!> m by n band matrix, with kl sub-diagonals and ku super-diagonals.
!> 
Parameters
[in]TRANS
!>          TRANS is CHARACTER*1
!>           On entry, TRANS specifies the operation to be performed as
!>           follows:
!>
!>              TRANS = 'N' or 'n'   y := alpha*A*x + beta*y.
!>
!>              TRANS = 'T' or 't'   y := alpha*A**T*x + beta*y.
!>
!>              TRANS = 'C' or 'c'   y := alpha*A**H*x + beta*y.
!> 
[in]M
!>          M is INTEGER
!>           On entry, M specifies the number of rows of the matrix A.
!>           M must be at least zero.
!> 
[in]N
!>          N is INTEGER
!>           On entry, N specifies the number of columns of the matrix A.
!>           N must be at least zero.
!> 
[in]KL
!>          KL is INTEGER
!>           On entry, KL specifies the number of sub-diagonals of the
!>           matrix A. KL must satisfy  0 .le. KL.
!> 
[in]KU
!>          KU is INTEGER
!>           On entry, KU specifies the number of super-diagonals of the
!>           matrix A. KU must satisfy  0 .le. KU.
!> 
[in]ALPHA
!>          ALPHA is COMPLEX
!>           On entry, ALPHA specifies the scalar alpha.
!> 
[in]A
!>          A is COMPLEX array, dimension ( LDA, N )
!>           Before entry, the leading ( kl + ku + 1 ) by n part of the
!>           array A must contain the matrix of coefficients, supplied
!>           column by column, with the leading diagonal of the matrix in
!>           row ( ku + 1 ) of the array, the first super-diagonal
!>           starting at position 2 in row ku, the first sub-diagonal
!>           starting at position 1 in row ( ku + 2 ), and so on.
!>           Elements in the array A that do not correspond to elements
!>           in the band matrix (such as the top left ku by ku triangle)
!>           are not referenced.
!>           The following program segment will transfer a band matrix
!>           from conventional full matrix storage to band storage:
!>
!>                 DO 20, J = 1, N
!>                    K = KU + 1 - J
!>                    DO 10, I = MAX( 1, J - KU ), MIN( M, J + KL )
!>                       A( K + I, J ) = matrix( I, J )
!>              10    CONTINUE
!>              20 CONTINUE
!> 
[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
!>           ( kl + ku + 1 ).
!> 
[in]X
!>          X is COMPLEX array, dimension at least
!>           ( 1 + ( n - 1 )*abs( INCX ) ) when TRANS = 'N' or 'n'
!>           and at least
!>           ( 1 + ( m - 1 )*abs( INCX ) ) otherwise.
!>           Before entry, the incremented array X must contain the
!>           vector x.
!> 
[in]INCX
!>          INCX is INTEGER
!>           On entry, INCX specifies the increment for the elements of
!>           X. INCX must not be zero.
!> 
[in]BETA
!>          BETA is COMPLEX
!>           On entry, BETA specifies the scalar beta. When BETA is
!>           supplied as zero then Y need not be set on input.
!> 
[in,out]Y
!>          Y is COMPLEX array, dimension at least
!>           ( 1 + ( m - 1 )*abs( INCY ) ) when TRANS = 'N' or 'n'
!>           and at least
!>           ( 1 + ( n - 1 )*abs( INCY ) ) otherwise.
!>           Before entry, the incremented array Y must contain the
!>           vector y. On exit, Y is overwritten by the updated vector y.
!>           If either m or n is zero, then Y not referenced and the function
!>           performs a quick return.
!> 
[in]INCY
!>          INCY is INTEGER
!>           On entry, INCY specifies the increment for the elements of
!>           Y. INCY 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 188 of file cgbmv.f.

190*
191* -- Reference BLAS level2 routine --
192* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
193* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
194*
195* .. Scalar Arguments ..
196 COMPLEX ALPHA,BETA
197 INTEGER INCX,INCY,KL,KU,LDA,M,N
198 CHARACTER TRANS
199* ..
200* .. Array Arguments ..
201 COMPLEX A(LDA,*),X(*),Y(*)
202* ..
203*
204* =====================================================================
205*
206* .. Parameters ..
207 COMPLEX ONE
208 parameter(one= (1.0e+0,0.0e+0))
209 COMPLEX ZERO
210 parameter(zero= (0.0e+0,0.0e+0))
211* ..
212* .. Local Scalars ..
213 COMPLEX TEMP
214 INTEGER I,INFO,IX,IY,J,JX,JY,K,KUP1,KX,KY,LENX,LENY
215 LOGICAL NOCONJ
216* ..
217* .. External Functions ..
218 LOGICAL LSAME
219 EXTERNAL lsame
220* ..
221* .. External Subroutines ..
222 EXTERNAL xerbla
223* ..
224* .. Intrinsic Functions ..
225 INTRINSIC conjg,max,min
226* ..
227*
228* Test the input parameters.
229*
230 info = 0
231 IF (.NOT.lsame(trans,'N') .AND. .NOT.lsame(trans,'T') .AND.
232 + .NOT.lsame(trans,'C')) THEN
233 info = 1
234 ELSE IF (m.LT.0) THEN
235 info = 2
236 ELSE IF (n.LT.0) THEN
237 info = 3
238 ELSE IF (kl.LT.0) THEN
239 info = 4
240 ELSE IF (ku.LT.0) THEN
241 info = 5
242 ELSE IF (lda.LT. (kl+ku+1)) THEN
243 info = 8
244 ELSE IF (incx.EQ.0) THEN
245 info = 10
246 ELSE IF (incy.EQ.0) THEN
247 info = 13
248 END IF
249 IF (info.NE.0) THEN
250 CALL xerbla('CGBMV ',info)
251 RETURN
252 END IF
253*
254* Quick return if possible.
255*
256 IF ((m.EQ.0) .OR. (n.EQ.0) .OR.
257 + ((alpha.EQ.zero).AND. (beta.EQ.one))) RETURN
258*
259 noconj = lsame(trans,'T')
260*
261* Set LENX and LENY, the lengths of the vectors x and y, and set
262* up the start points in X and Y.
263*
264 IF (lsame(trans,'N')) THEN
265 lenx = n
266 leny = m
267 ELSE
268 lenx = m
269 leny = n
270 END IF
271 IF (incx.GT.0) THEN
272 kx = 1
273 ELSE
274 kx = 1 - (lenx-1)*incx
275 END IF
276 IF (incy.GT.0) THEN
277 ky = 1
278 ELSE
279 ky = 1 - (leny-1)*incy
280 END IF
281*
282* Start the operations. In this version the elements of A are
283* accessed sequentially with one pass through the band part of A.
284*
285* First form y := beta*y.
286*
287 IF (beta.NE.one) THEN
288 IF (incy.EQ.1) THEN
289 IF (beta.EQ.zero) THEN
290 DO 10 i = 1,leny
291 y(i) = zero
292 10 CONTINUE
293 ELSE
294 DO 20 i = 1,leny
295 y(i) = beta*y(i)
296 20 CONTINUE
297 END IF
298 ELSE
299 iy = ky
300 IF (beta.EQ.zero) THEN
301 DO 30 i = 1,leny
302 y(iy) = zero
303 iy = iy + incy
304 30 CONTINUE
305 ELSE
306 DO 40 i = 1,leny
307 y(iy) = beta*y(iy)
308 iy = iy + incy
309 40 CONTINUE
310 END IF
311 END IF
312 END IF
313 IF (alpha.EQ.zero) RETURN
314 kup1 = ku + 1
315 IF (lsame(trans,'N')) THEN
316*
317* Form y := alpha*A*x + y.
318*
319 jx = kx
320 IF (incy.EQ.1) THEN
321 DO 60 j = 1,n
322 temp = alpha*x(jx)
323 k = kup1 - j
324 DO 50 i = max(1,j-ku),min(m,j+kl)
325 y(i) = y(i) + temp*a(k+i,j)
326 50 CONTINUE
327 jx = jx + incx
328 60 CONTINUE
329 ELSE
330 DO 80 j = 1,n
331 temp = alpha*x(jx)
332 iy = ky
333 k = kup1 - j
334 DO 70 i = max(1,j-ku),min(m,j+kl)
335 y(iy) = y(iy) + temp*a(k+i,j)
336 iy = iy + incy
337 70 CONTINUE
338 jx = jx + incx
339 IF (j.GT.ku) ky = ky + incy
340 80 CONTINUE
341 END IF
342 ELSE
343*
344* Form y := alpha*A**T*x + y or y := alpha*A**H*x + y.
345*
346 jy = ky
347 IF (incx.EQ.1) THEN
348 DO 110 j = 1,n
349 temp = zero
350 k = kup1 - j
351 IF (noconj) THEN
352 DO 90 i = max(1,j-ku),min(m,j+kl)
353 temp = temp + a(k+i,j)*x(i)
354 90 CONTINUE
355 ELSE
356 DO 100 i = max(1,j-ku),min(m,j+kl)
357 temp = temp + conjg(a(k+i,j))*x(i)
358 100 CONTINUE
359 END IF
360 y(jy) = y(jy) + alpha*temp
361 jy = jy + incy
362 110 CONTINUE
363 ELSE
364 DO 140 j = 1,n
365 temp = zero
366 ix = kx
367 k = kup1 - j
368 IF (noconj) THEN
369 DO 120 i = max(1,j-ku),min(m,j+kl)
370 temp = temp + a(k+i,j)*x(ix)
371 ix = ix + incx
372 120 CONTINUE
373 ELSE
374 DO 130 i = max(1,j-ku),min(m,j+kl)
375 temp = temp + conjg(a(k+i,j))*x(ix)
376 ix = ix + incx
377 130 CONTINUE
378 END IF
379 y(jy) = y(jy) + alpha*temp
380 jy = jy + incy
381 IF (j.GT.ku) kx = kx + incx
382 140 CONTINUE
383 END IF
384 END IF
385*
386 RETURN
387*
388* End of CGBMV
389*
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: