ZGEMV performs one of the matrix-vector operations y := alphaAx + betay, or y := alphaATx + betay, or y := alpha*AHx + betay, where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
character(len=1), | intent(in) | :: | trans | |||
integer(kind=ilp), | intent(in) | :: | m | |||
integer(kind=ilp), | intent(in) | :: | n | |||
complex(kind=dp), | intent(in) | :: | alpha | |||
complex(kind=dp), | intent(in) | :: | a(lda,*) | |||
integer(kind=ilp), | intent(in) | :: | lda | |||
complex(kind=dp), | intent(in) | :: | x(*) | |||
integer(kind=ilp), | intent(in) | :: | incx | |||
complex(kind=dp), | intent(in) | :: | beta | |||
complex(kind=dp), | intent(inout) | :: | y(*) | |||
integer(kind=ilp), | intent(in) | :: | incy |
pure subroutine stdlib_zgemv(trans,m,n,alpha,a,lda,x,incx,beta,y,incy) !! ZGEMV 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 matrix. ! -- reference blas level2 routine -- ! -- reference blas is a software package provided by univ. of tennessee, -- ! -- univ. of california berkeley, univ. of colorado denver and nag ltd..-- ! Scalar Arguments complex(dp), intent(in) :: alpha, beta integer(ilp), intent(in) :: incx, incy, lda, m, n character, intent(in) :: trans ! Array Arguments complex(dp), intent(in) :: a(lda,*), x(*) complex(dp), intent(inout) :: y(*) ! ===================================================================== ! Local Scalars complex(dp) :: temp integer(ilp) :: i, info, ix, iy, j, jx, jy, kx, ky, lenx, leny logical(lk) :: noconj ! Intrinsic Functions intrinsic :: conjg,max ! test the input parameters. info = 0 if (.not.stdlib_lsame(trans,'N') .and. .not.stdlib_lsame(trans,'T') & .and..not.stdlib_lsame(trans,'C')) then info = 1 else if (m<0) then info = 2 else if (n<0) then info = 3 else if (lda<max(1,m)) then info = 6 else if (incx==0) then info = 8 else if (incy==0) then info = 11 end if if (info/=0) then call stdlib_xerbla('ZGEMV ',info) return end if ! quick return if possible. if ((m==0) .or. (n==0) .or.((alpha==czero).and. (beta==cone))) return noconj = stdlib_lsame(trans,'T') ! set lenx and leny, the lengths of the vectors x and y, and set ! up the start points in x and y. if (stdlib_lsame(trans,'N')) then lenx = n leny = m else lenx = m leny = n end if if (incx>0) then kx = 1 else kx = 1 - (lenx-1)*incx end if if (incy>0) then ky = 1 else ky = 1 - (leny-1)*incy end if ! start the operations. in this version the elements of a are ! accessed sequentially with cone pass through a. ! first form y := beta*y. if (beta/=cone) then if (incy==1) then if (beta==czero) then do i = 1,leny y(i) = czero end do else do i = 1,leny y(i) = beta*y(i) end do end if else iy = ky if (beta==czero) then do i = 1,leny y(iy) = czero iy = iy + incy end do else do i = 1,leny y(iy) = beta*y(iy) iy = iy + incy end do end if end if end if if (alpha==czero) return if (stdlib_lsame(trans,'N')) then ! form y := alpha*a*x + y. jx = kx if (incy==1) then do j = 1,n temp = alpha*x(jx) do i = 1,m y(i) = y(i) + temp*a(i,j) end do jx = jx + incx end do 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 end do jx = jx + incx end do end if else ! form y := alpha*a**t*x + y or y := alpha*a**h*x + y. jy = ky if (incx==1) then do j = 1,n temp = czero if (noconj) then do i = 1,m temp = temp + a(i,j)*x(i) end do else do i = 1,m temp = temp + conjg(a(i,j))*x(i) end do end if y(jy) = y(jy) + alpha*temp jy = jy + incy end do else do j = 1,n temp = czero ix = kx if (noconj) then do i = 1,m temp = temp + a(i,j)*x(ix) ix = ix + incx end do else do i = 1,m temp = temp + conjg(a(i,j))*x(ix) ix = ix + incx end do end if y(jy) = y(jy) + alpha*temp jy = jy + incy end do end if end if return end subroutine stdlib_zgemv