ddot(3f) - [BLAS:DOUBLE_BLAS_LEVEL1] forms the dot product of two vectors.
double precision function ddot(n,dx,incx,dy,incy)
.. Scalar Arguments ..
integer,intent(in) :: incx,incy,n
..
.. Array Arguments ..
double precision,intent(in) :: dx(*),dy(*)
..
DDOT forms the dot product of two vectors.
uses unrolled loops for increments equal to one.
N
N is INTEGER
number of elements in input vector(s)
DX
DX is DOUBLE PRECISION array, dimension ( 1 + ( N - 1 )*abs( INCX ) )
INCX
INCX is INTEGER
storage spacing between elements of DX
DY
DY is DOUBLE PRECISION array, dimension ( 1 + ( N - 1 )*abs( INCY ) )
INCY
INCY is INTEGER
storage spacing between elements of DY
date:November 2017
FURTHER DETAILS
jack dongarra, linpack, 3/11/78.
modified 12/3/93, array(1) declarations changed to array(*)
Online html documentation available at
http://www.netlib.org/lapack/explore-html/
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer, | intent(in) | :: | n | |||
double precision, | intent(in) | :: | dx(*) | |||
integer, | intent(in) | :: | incx | |||
double precision, | intent(in) | :: | dy(*) | |||
integer, | intent(in) | :: | incy |
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
double precision, | public | :: | dtemp | ||||
integer, | public | :: | i | ||||
integer, | public | :: | ix | ||||
integer, | public | :: | iy | ||||
integer, | public | :: | m | ||||
integer, | public | :: | mp1 |
pure double precision function ddot(n,dx,incx,dy,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 ..
integer,intent(in) :: incx,incy,n
! ..
! .. Array Arguments ..
double precision,intent(in) :: dx(*),dy(*)
! ..
!
! =====================================================================
!
! .. Local Scalars ..
double precision dtemp
integer i,ix,iy,m,mp1
! ..
! .. Intrinsic Functions ..
intrinsic mod
! ..
ddot = 0.0d0
dtemp = 0.0d0
if (n.le.0) return
if (incx.eq.1 .and. incy.eq.1) then
!
! code for both increments equal to 1
!
!
! clean-up loop
!
m = mod(n,5)
if (m.ne.0) then
do i = 1,m
dtemp = dtemp + dx(i)*dy(i)
enddo
if (n.lt.5) then
ddot=dtemp
return
endif
endif
mp1 = m + 1
do i = mp1,n,5
dtemp = dtemp + dx(i)*dy(i) + dx(i+1)*dy(i+1) + dx(i+2)*dy(i+2) + dx(i+3)*dy(i+3) + dx(i+4)*dy(i+4)
enddo
else
!
! code for unequal increments or equal increments
! not equal to 1
!
ix = 1
iy = 1
if (incx.lt.0) ix = (-n+1)*incx + 1
if (incy.lt.0) iy = (-n+1)*incy + 1
do i = 1,n
dtemp = dtemp + dx(ix)*dy(iy)
ix = ix + incx
iy = iy + incy
enddo
endif
ddot = dtemp
end function ddot