sdsdot(3f) - [BLAS:SINGLE_BLAS_LEVEL1] Compute the inner
product of two vectors with extended precision accumulation.
SDSDOT := SUM SX * SY (accumulated double precision, returned single)
real function sdsdot(n,sb,sx,incx,sy,incy)
.. Scalar Arguments ..
real,intent(in) :: sb
integer,intent(in) :: incx,incy,n
..
.. Array Arguments ..
real,intent(in) :: sx(*),sy(*)
..
Compute the inner product of two vectors with extended
precision accumulation.
Returns S.P. result with dot product accumulated in D.P.
SDSDOT = SB + sum for I = 0 to N-1 of SX(LX+I*INCX)*SY(LY+I*INCY),
where LX = 1 if INCX .GE. 0, else LX = 1+(1-N)*INCX, and LY is
defined in a similar way using INCY.
N
N is INTEGER
number of elements in input vector(s)
SB
SB is REAL
single precision scalar to be added to inner product
SX
SX is REAL array, dimension ( 1 + ( N - 1 )*abs( INCX ) )
single precision vector with N elements
INCX
INCX is INTEGER
storage spacing between elements of SX
SY
SY is REAL array, dimension ( 1 + ( N - 1 )*abs( INCX ) )
single precision vector with N elements
INCY
INCY is INTEGER
storage spacing between elements of SY
Kincaid, D. R., (U. of Texas), Krogh, F. T., (JPL)
Univ. of Tennessee
date:November 2017
FURTHER DETAILS
REFERENCES
C. L. Lawson, R. J. Hanson, D. R. Kincaid and F. T.
Krogh, Basic linear algebra subprograms for Fortran
usage, Algorithm No. 539, Transactions on Mathematical
Software 5, 3 (September 1979), pp. 308-323.
REVISION HISTORY (YYMMDD)
791001 DATE WRITTEN
890531 Changed all specific intrinsics to generic. (WRB)
890831 Modified array declarations. (WRB)
890831 REVISION DATE from Version 3.2
891214 Prologue converted to Version 4.0 format. (BAB)
920310 Corrected definition of LX in DESCRIPTION. (WRB)
920501 Reformatted the REFERENCES section. (WRB)
070118 Reformat to LAPACK coding style
Online html documentation available at
http://www.netlib.org/lapack/explore-html/
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | |||
real, | intent(in) | :: | sb | |||
real, | intent(in) | :: | sx(*) | |||
integer, | intent(in) | :: | incx | |||
real, | intent(in) | :: | sy(*) | |||
integer, | intent(in) | :: | incy |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
double precision, | public | :: | dsdot | ||||
integer, | public | :: | i | ||||
integer, | public | :: | kx | ||||
integer, | public | :: | ky | ||||
integer, | public | :: | ns |
pure real function sdsdot(n,sb,sx,incx,sy,incy)
implicit none
!
! -- Reference BLAS level1 routine (version 3.8.0) --
! -- Reference BLAS is a software package provided by Univ. of Tennessee, --
! -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
! November 2017
!
! .. Scalar Arguments ..
real,intent(in) :: sb
integer,intent(in) :: incx,incy,n
! ..
! .. Array Arguments ..
real,intent(in) :: sx(*),sy(*)
! ..
! .. Local Scalars ..
double precision dsdot
integer i,kx,ky,ns
! ..
! .. Intrinsic Functions ..
intrinsic dble
! ..
dsdot = sb
if (n.le.0) then
sdsdot = real(dsdot)
return
endif
if (incx.eq.incy .and. incx.gt.0) then
!
! Code for equal and positive increments.
!
ns = n*incx
do i = 1,ns,incx
dsdot = dsdot + dble(sx(i))*dble(sy(i))
enddo
else
!
! Code for unequal or nonpositive increments.
!
kx = 1
ky = 1
if (incx.lt.0) kx = 1 + (1-n)*incx
if (incy.lt.0) ky = 1 + (1-n)*incy
do i = 1,n
dsdot = dsdot + dble(sx(kx))*dble(sy(ky))
kx = kx + incx
ky = ky + incy
enddo
endif
sdsdot = real(dsdot)
end function sdsdot