SSPMV performs the matrix-vector operation y := alphaAx + beta*y, where alpha and beta are scalars, x and y are n element vectors and A is an n by n symmetric matrix, supplied in packed form.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
character(len=1), | intent(in) | :: | uplo | |||
integer(kind=ilp), | intent(in) | :: | n | |||
real(kind=sp), | intent(in) | :: | alpha | |||
real(kind=sp), | intent(in) | :: | ap(*) | |||
real(kind=sp), | intent(in) | :: | x(*) | |||
integer(kind=ilp), | intent(in) | :: | incx | |||
real(kind=sp), | intent(in) | :: | beta | |||
real(kind=sp), | intent(inout) | :: | y(*) | |||
integer(kind=ilp), | intent(in) | :: | incy |
pure subroutine stdlib_sspmv(uplo,n,alpha,ap,x,incx,beta,y,incy) !! SSPMV performs the matrix-vector operation !! y := alpha*A*x + beta*y, !! where alpha and beta are scalars, x and y are n element vectors and !! A is an n by n symmetric matrix, supplied in packed form. ! -- 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 real(sp), intent(in) :: alpha, beta integer(ilp), intent(in) :: incx, incy, n character, intent(in) :: uplo ! Array Arguments real(sp), intent(in) :: ap(*), x(*) real(sp), intent(inout) :: y(*) ! ===================================================================== ! Local Scalars real(sp) :: temp1, temp2 integer(ilp) :: i, info, ix, iy, j, jx, jy, k, kk, kx, ky ! test the input parameters. info = 0 if (.not.stdlib_lsame(uplo,'U') .and. .not.stdlib_lsame(uplo,'L')) then info = 1 else if (n<0) then info = 2 else if (incx==0) then info = 6 else if (incy==0) then info = 9 end if if (info/=0) then call stdlib_xerbla('SSPMV ',info) return end if ! quick return if possible. if ((n==0) .or. ((alpha==zero).and. (beta==one))) return ! set up the start points in x and y. if (incx>0) then kx = 1 else kx = 1 - (n-1)*incx end if if (incy>0) then ky = 1 else ky = 1 - (n-1)*incy end if ! start the operations. in this version the elements of the array ap ! are accessed sequentially with one pass through ap. ! first form y := beta*y. if (beta/=one) then if (incy==1) then if (beta==zero) then do i = 1,n y(i) = zero end do else do i = 1,n y(i) = beta*y(i) end do end if else iy = ky if (beta==zero) then do i = 1,n y(iy) = zero iy = iy + incy end do else do i = 1,n y(iy) = beta*y(iy) iy = iy + incy end do end if end if end if if (alpha==zero) return kk = 1 if (stdlib_lsame(uplo,'U')) then ! form y when ap contains the upper triangle. if ((incx==1) .and. (incy==1)) then do j = 1,n temp1 = alpha*x(j) temp2 = zero k = kk do i = 1,j - 1 y(i) = y(i) + temp1*ap(k) temp2 = temp2 + ap(k)*x(i) k = k + 1 end do y(j) = y(j) + temp1*ap(kk+j-1) + alpha*temp2 kk = kk + j end do else jx = kx jy = ky do j = 1,n temp1 = alpha*x(jx) temp2 = zero ix = kx iy = ky do k = kk,kk + j - 2 y(iy) = y(iy) + temp1*ap(k) temp2 = temp2 + ap(k)*x(ix) ix = ix + incx iy = iy + incy end do y(jy) = y(jy) + temp1*ap(kk+j-1) + alpha*temp2 jx = jx + incx jy = jy + incy kk = kk + j end do end if else ! form y when ap contains the lower triangle. if ((incx==1) .and. (incy==1)) then do j = 1,n temp1 = alpha*x(j) temp2 = zero y(j) = y(j) + temp1*ap(kk) k = kk + 1 do i = j + 1,n y(i) = y(i) + temp1*ap(k) temp2 = temp2 + ap(k)*x(i) k = k + 1 end do y(j) = y(j) + alpha*temp2 kk = kk + (n-j+1) end do else jx = kx jy = ky do j = 1,n temp1 = alpha*x(jx) temp2 = zero y(jy) = y(jy) + temp1*ap(kk) ix = jx iy = jy do k = kk + 1,kk + n - j ix = ix + incx iy = iy + incy y(iy) = y(iy) + temp1*ap(k) temp2 = temp2 + ap(k)*x(ix) end do y(jy) = y(jy) + alpha*temp2 jx = jx + incx jy = jy + incy kk = kk + (n-j+1) end do end if end if return end subroutine stdlib_sspmv