Section: BLAS/LAPACK - DDOT Timing#
Adapted from: gjbex/Fortran-MOOC
This program demonstrates copying vectors using BLAS/LAPACK in Fortran.#
program ddot_timing
use, intrinsic :: iso_fortran_env, only : error_unit, DP => REAL64
implicit none
integer :: n, runs, run, i, istat
real(kind=DP), dimension(:), allocatable :: u, v
real(kind=DP) :: scale_u, scale_v, r, r_ddot
real(kind=DP) :: ddot_time, loop_time
real :: start_time, end_time
interface
real(kind=DP) function ddot(n, dx, incx, dy, incy)
import :: DP
integer :: n, incx, incy
real(kind=DP), dimension(*) :: dx, dy
end function ddot
end interface
call get_arguments(runs, n)
allocate (u(n), v(n), stat=istat)
if (istat /= 0) then
write (unit=error_unit, fmt='(A)') 'error: can not allocate vectors'
stop 2
end if
do run = 1, runs
print *, "Received command line arguments..."
print '(A, I10)', 'Runs: ', runs
print '(A, I10)', 'Array size: ', n
print *, 'Filling arrays u, v with data...'
scale_u = 1.17*run
scale_v = 0.17*run
do i = 1, n
u(i) = scale_u*sqrt(real(i))/n
v(i) = scale_v*sqrt(real(i))/n
end do
print *, 'Starting DDOT run...'
call cpu_time(start_time)
r_ddot = ddot(n, u, 1, v, 1)
call cpu_time(end_time)
ddot_time = end_time - start_time
call cpu_time(start_time)
print *, 'Starting Loop run...'
r = 0.0_DP
do i = 1, n
r = r + u(i)*v(i)
end do
call cpu_time(end_time)
loop_time = end_time - start_time
print '(A, F10.6, A, F10.6)', 'ddot time: ', ddot_time, &
', loop time: ', loop_time
end do
deallocate (u, v)
contains
subroutine get_arguments(runs, n)
use, intrinsic :: iso_fortran_env, only : error_unit
implicit none
integer, intent(out) :: runs, n
integer :: istat
character(len=1024) :: buffer, msg
runs = 1
if (command_argument_count() >= 1) then
call get_command_argument(1, buffer)
read (buffer, fmt=*, iostat=istat, iomsg=msg) runs
if (istat /= 0) then
write (unit=error_unit, fmt='(2A)') 'error: ', trim(msg)
stop 1
end if
end if
n = 50000000
if (command_argument_count() >= 2) then
call get_command_argument(2, buffer)
read (buffer, fmt=*) n
if (istat /= 0) then
write (unit=error_unit, fmt='(2A)') 'error: ', trim(msg)
stop 1
end if
end if
end subroutine get_arguments
end program ddot_timing
The above program is compiled and run using Fortran Package Manager (fpm):
Build the Program using FPM (Fortran Package Manager)#
import os
root_dir = ""
root_dir = os.getcwd()
Since the code makes use of the LAPACK library, the following FPM configuration file (fpm.toml) was used:
name = "Section_BLAS_LAPACK_DDOT_Timing"
[build]
auto-executables = true
auto-tests = true
auto-examples = true
link = ["blas", "lapack"]
[install]
library = false
[[executable]]
name="Section_BLAS_LAPACK_DDOT_Timing"
source-dir="app"
main="section_blas_lapack_ddot_timing.f90"
code_dir = root_dir + "/" + "Fortran_Code/Section_BLAS_LAPACK_DDOT_Timing"
os.chdir(code_dir)
build_status = os.system("fpm build 2>/dev/null")
Run the Program using FPM (Fortran Package Manager)#
exec_status = \
os.system("fpm run 2>/dev/null -- 1 1000")
Received command line arguments...
Runs: 1
Array size: 1000
Filling arrays u, v with data...
Starting DDOT run...
Starting Loop run...
ddot time: 0.000030, loop time: 0.000006
exec_status = \
os.system("fpm run 2>/dev/null -- 1 10000")
Received command line arguments...
Runs: 1
Array size: 10000
Filling arrays u, v with data...
Starting DDOT run...
Starting Loop run...
ddot time: 0.000033, loop time: 0.000022
exec_status = \
os.system("fpm run 2>/dev/null -- 1 100000")
Received command line arguments...
Runs: 1
Array size: 100000
Filling arrays u, v with data...
Starting DDOT run...
Starting Loop run...
ddot time: 0.000224, loop time: 0.000427
exec_status = \
os.system("fpm run 2>/dev/null -- 1 1000000")
Received command line arguments...
Runs: 1
Array size: 1000000
Filling arrays u, v with data...
Starting DDOT run...
Starting Loop run...
ddot time: 0.002128, loop time: 0.002818
exec_status = \
os.system("fpm run 2>/dev/null -- 1 10000000")
Received command line arguments...
Runs: 1
Array size: 10000000
Filling arrays u, v with data...
Starting DDOT run...
Starting Loop run...
ddot time: 0.017535, loop time: 0.031079
exec_status = \
os.system("fpm run 2>/dev/null -- 1 100000000")
Received command line arguments...
Runs: 1
Array size: 100000000
Filling arrays u, v with data...
Starting DDOT run...
Starting Loop run...
ddot time: 0.244720, loop time: 0.387014
exec_status = \
os.system("fpm run 2>/dev/null")
Received command line arguments...
Runs: 1
Array size: 50000000
Filling arrays u, v with data...
Starting DDOT run...
Starting Loop run...
ddot time: 0.100836, loop time: 0.177330