dgemv(3f) - [BLAS:DOUBLE_BLAS_LEVEL2]
subroutine dgemv(trans,m,n,alpha,a,lda,x,incx,beta,y,incy)
.. Scalar Arguments ..
double precision,intent(in) :: alpha,beta
integer,intent(in) :: incx,incy,lda,m,n
character,intent(in) :: trans
..
.. Array Arguments ..
double precision,intent(in) :: a(lda,*),x(*)
double precision,intent(inout) :: y(*)
..
DGEMV performs one of the matrix-vector operations
y := alpha*A*x + beta*y, or y := alpha*A**T*x + beta*y,
where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.
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**T*x + beta*y.
M
M is INTEGER
On entry, M specifies the number of rows of the matrix A.
M must be at least zero.
N
N is INTEGER
On entry, N specifies the number of columns of the matrix A.
N must be at least zero.
ALPHA
ALPHA is DOUBLE PRECISION.
On entry, ALPHA specifies the scalar alpha.
A
A is DOUBLE PRECISION array, dimension ( LDA, N )
Before entry, the leading m by n part of the array A must
contain the matrix of coefficients.
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
max( 1, m ).
X
X is DOUBLE PRECISION 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.
INCX
INCX is INTEGER
On entry, INCX specifies the increment for the elements of
X. INCX must not be zero.
BETA
BETA is DOUBLE PRECISION.
On entry, BETA specifies the scalar beta. When BETA is
supplied as zero then Y need not be set on input.
Y
Y is DOUBLE PRECISION 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 with BETA non-zero, the incremented array Y
must contain the vector y. On exit, Y is overwritten by the
updated vector y.
INCY
INCY is INTEGER
On entry, INCY specifies the increment for the elements of
Y. INCY must not be zero.
date:December 2016
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.
Online html documentation available at
http://www.netlib.org/lapack/explore-html/
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
character(len=1), | intent(in) | :: | trans | |||
integer, | intent(in) | :: | m | |||
integer, | intent(in) | :: | n | |||
double precision, | intent(in) | :: | alpha | |||
double precision, | intent(in) | :: | a(lda,*) | |||
integer, | intent(in) | :: | lda | |||
double precision, | intent(in) | :: | x(*) | |||
integer, | intent(in) | :: | incx | |||
double precision, | intent(in) | :: | beta | |||
double precision, | intent(inout) | :: | y(*) | |||
integer, | intent(in) | :: | incy |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public | :: | i | ||||
integer, | public | :: | info | ||||
integer, | public | :: | ix | ||||
integer, | public | :: | iy | ||||
integer, | public | :: | j | ||||
integer, | public | :: | jx | ||||
integer, | public | :: | jy | ||||
integer, | public | :: | kx | ||||
integer, | public | :: | ky | ||||
integer, | public | :: | lenx | ||||
integer, | public | :: | leny | ||||
double precision, | public, | parameter | :: | one | = | 1.0d+0 | |
double precision, | public | :: | temp | ||||
double precision, | public, | parameter | :: | zero | = | 0.0d+0 |
subroutine dgemv(trans,m,n,alpha,a,lda,x,incx,beta,y,incy)
implicit none
!
! -- Reference BLAS level2 routine (version 3.7.0) --
! -- Reference BLAS is a software package provided by Univ. of Tennessee, --
! -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
! December 2016
!
! .. Scalar Arguments ..
double precision,intent(in) :: alpha,beta
integer,intent(in) :: incx,incy,lda,m,n
character,intent(in) :: trans
! ..
! .. Array Arguments ..
double precision,intent(in) :: a(lda,*),x(*)
double precision,intent(inout) :: y(*)
! ..
!
! =====================================================================
!
! .. Parameters ..
double precision one,zero
parameter (one=1.0d+0,zero=0.0d+0)
! ..
! .. Local Scalars ..
double precision temp
integer i,info,ix,iy,j,jx,jy,kx,ky,lenx,leny
! ..
! .. External Functions ..
! ..
! .. External Subroutines ..
! ..
! .. Intrinsic Functions ..
intrinsic max
! ..
!
! Test the input parameters.
!
info = 0
if (.not.lsame(trans,'N') .and. .not.lsame(trans,'T') .and. .not.lsame(trans,'C')) then
info = 1
else if (m.lt.0) then
info = 2
else if (n.lt.0) then
info = 3
else if (lda.lt.max(1,m)) then
info = 6
else if (incx.eq.0) then
info = 8
else if (incy.eq.0) then
info = 11
endif
if (info.ne.0) then
call xerbla('DGEMV ',info)
return
endif
!
! Quick return if possible.
!
if ((m.eq.0) .or. (n.eq.0) .or. ((alpha.eq.zero).and. (beta.eq.one))) return
!
! Set LENX and LENY, the lengths of the vectors x and y, and set
! up the start points in X and Y.
!
if (lsame(trans,'N')) then
lenx = n
leny = m
else
lenx = m
leny = n
endif
if (incx.gt.0) then
kx = 1
else
kx = 1 - (lenx-1)*incx
endif
if (incy.gt.0) then
ky = 1
else
ky = 1 - (leny-1)*incy
endif
!
! Start the operations. In this version the elements of A are
! accessed sequentially with one pass through A.
!
! First form y := beta*y.
!
if (beta.ne.one) then
if (incy.eq.1) then
if (beta.eq.zero) then
y(1:leny) = zero
else
y(1:leny) = beta*y(1:leny)
endif
else
iy = ky
if (beta.eq.zero) then
do i = 1,leny
y(iy) = zero
iy = iy + incy
enddo
else
do i = 1,leny
y(iy) = beta*y(iy)
iy = iy + incy
enddo
endif
endif
endif
if (alpha.eq.zero) return
if (lsame(trans,'N')) then
!
! Form y := alpha*A*x + y.
!
jx = kx
if (incy.eq.1) then
do j = 1,n
temp = alpha*x(jx)
y(1:m) = y(1:m) + temp*a(1:m,j)
jx = jx + incx
enddo
else
do j = 1,n
temp = alpha*x(jx)
iy = ky
do i = 1,m
y(iy) = y(iy) + temp*a(i,j)
iy = iy + incy
enddo
jx = jx + incx
enddo
endif
else
!
! Form y := alpha*A**T*x + y.
!
jy = ky
if (incx.eq.1) then
do j = 1,n
temp = zero
do i = 1,m
temp = temp + a(i,j)*x(i)
enddo
y(jy) = y(jy) + alpha*temp
jy = jy + incy
enddo
else
do j = 1,n
temp = zero
ix = kx
do i = 1,m
temp = temp + a(i,j)*x(ix)
ix = ix + incx
enddo
y(jy) = y(jy) + alpha*temp
jy = jy + incy
enddo
endif
endif
!
end subroutine dgemv