You cannot select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
psblas3/test/computational_routines/gedot/gedot.f90

227 lines
8.1 KiB
Fortran

program main
use psb_gedot_test
use psb_base_mod
implicit none
! MPI variables
integer(psb_ipk_) :: my_rank, np
! Communicator variable
type(psb_ctxt_type) :: ctxt
! parameters array
character(len=64) :: x(4),y(4)
real(psb_dpk_) :: alpha(3), beta(3)
integer(psb_ipk_) :: arr_size
integer(psb_ipk_) :: tests_number, count
! cycle indexes variables
integer(psb_ipk_) :: i,j,k,h
integer(psb_ipk_) :: info, ret, unit
! time stats variables
character(len=8) :: date ! YYYYMMDD
character(len=10) :: time ! HHMMSS.sss
character(len=5) :: zones ! Time zone
integer :: values(8)
! others
character(len=:), allocatable :: output_file_name
integer(psb_ipk_) :: last_percent_sp, last_percent_dp
! Initialize parameters
x(1) = "vectors/x1.mtx"
x(2) = "vectors/x2.mtx"
x(3) = "vectors/x3.mtx"
x(4) = "vectors/x4.mtx"
y(1) = "vectors/y1.mtx"
y(2) = "vectors/y2.mtx"
y(3) = "vectors/y3.mtx"
y(4) = "vectors/y4.mtx"
arr_size = 10000
tests_number = size(x) * size(y)
count = 0
last_percent_sp = -1
last_percent_dp = -1
call psb_init(ctxt)
call psb_info(ctxt,my_rank,np)
if(my_rank == psb_root_) then
! Setup logger output
if(np == 1) then
open(newunit=unit, file='psblas_gedot_test.log', status='replace', action='write', iostat=info)
else
open(newunit=unit, file='psblas_gedot_test.log', status='old', action='write', position='append', iostat=info)
end if
if (info /= 0) then
print *, 'Error opening output file.'
print *, "I/O Status Code:", info
stop
end if
psb_out_unit = unit
write(psb_out_unit,'(A,A)') 'Welcome to PSBLAS version: ',psb_version_string_
write(psb_out_unit,'(A)') 'This is the psb_gedot_test sample program'
write(psb_out_unit,'(A,I0)') 'Number of processes used in this computation: ', np
write(psb_out_unit,'(A)') ''
call generate_vectors(arr_size)
end if
call psb_bcast(ctxt,psb_out_unit)
call psb_barrier(ctxt)
if(my_rank == psb_root_) write(*,'(A)') "[INFO] Starting single precision computation..."
do i=1,size(x)
do j=1,size(y)
call psb_gedot_kernel(x_file=x(i), y_file=y(j), arr_size = arr_size, ctxt = ctxt, &
& ret = ret, output_file_name = output_file_name)
if(my_rank == psb_root_) then
count = count + 1
call print_progress(count, tests_number, last_percent_sp, "single precision")
call date_and_time(date, time, zones, values)
if(ret /= -1) then
! Success formatted output
write(psb_out_unit,'("[", I4.4,"-",I2.2,"-",I2.2," ",I2.2,":",I2.2,":",I2.2,"] ",&
& A,A,A,I0,A,I0,T110,A)') &
& values(1), values(2), values(3), values(5), values(6), values(7), &
& "Generation gedot single precision result file ", &
& output_file_name , ' ', count , "/", tests_number, "[OK]"
else
! Fail formatted output
write(psb_out_unit,'("[", I4.4,"-",I2.2,"-",I2.2," ",I2.2,":",I2.2,":",I2.2,"] ",&
& A,A,A,I0,A,I0,T110,A)') &
& values(1), values(2), values(3), values(5), values(6), values(7), &
& "Generation gedot single precision result file ", &
& output_file_name , ' ', count , "/", tests_number, "[FAIL]"
goto 9998
end if
end if
call psb_barrier(ctxt)
end do
end do
if(my_rank == psb_root_) write(*,'(A)') "[INFO] Single precision computation completed succesfully!"
if(my_rank == psb_root_) then
write(psb_out_unit, *) ''
count = 0
end if
if(my_rank == psb_root_) write(*,'(A)') "[INFO] Starting double precision check..."
call psb_barrier(ctxt)
! Here double precision comparison should be done
do i=1,size(x)
do j=1,size(y)
call psb_gedot_check(x_file=x(i), y_file=y(j), &
& arr_size = arr_size, ctxt = ctxt, ret = ret, output_file_name = output_file_name)
if(my_rank == psb_root_) then
count = count + 1
call print_progress(count, tests_number, last_percent_dp, "double precision")
call date_and_time(date, time, zones, values)
if(ret == 0) then
! Success formatted output
write(psb_out_unit,'("[", I4.4,"-",I2.2,"-",I2.2," ",I2.2,":",I2.2,":",I2.2,"] ",&
& A,A,A,I0,A,I0,T110,A)') &
& values(1), values(2), values(3), values(5), values(6), values(7), &
& "Double precision check on file ", &
& output_file_name , ' ', count , "/", tests_number, "[OK]"
else
! Fail formatted output
write(psb_out_unit,'("[", I4.4,"-",I2.2,"-",I2.2," ",I2.2,":",I2.2,":",I2.2,"] ",&
& A,A,A,I0,A,I0,T110,A)') &
& values(1), values(2), values(3), values(5), values(6), values(7), &
& "Double precision check on file ", &
& output_file_name , ' ', count , "/", tests_number, "[FAIL]"
write(psb_out_unit,'(A,I0)') "[ERROR] Error at element ", abs(ret)
goto 9999
end if
end if
call psb_barrier(ctxt)
end do
end do
if(my_rank == psb_root_) then
write(*,'(A)') "[INFO] Duble precision check completed succesfully!"
close(unit)
end if
call psb_exit(ctxt)
return
9998 continue
if(my_rank == psb_root_) then
close(unit)
write(*,'(A,I0,A,I0,A)') "[ERROR] Error in gedot single precision computation ", &
& count, "/", tests_number, " see log file for details"
end if
9999 continue
if(my_rank == psb_root_) then
close(unit)
write(*,'(A,I0,A,I0,A)') "[ERROR] Error in gedot double precision check ", &
& count, "/", tests_number, " see log file for details"
end if
call psb_exit(ctxt)
return
contains
subroutine print_progress(current, total, last_percent, label)
use iso_fortran_env, only: error_unit
implicit none
integer(psb_ipk_), intent(in) :: current, total
integer(psb_ipk_), intent(inout) :: last_percent
character(len=*), intent(in) :: label
integer(psb_ipk_) :: percent, filled, width, i
character(len=160) :: line
if (total <= 0) return
percent = int(real(current) / real(total) * 100.0)
if (percent == last_percent) return
last_percent = percent
width = 30
filled = int(real(percent) / 100.0 * width)
line = "[INFO] Progress " // trim(label) // ": ["
do i = 1, filled
line = trim(line) // "#"
end do
do i = filled + 1, width
line = trim(line) // "-"
end do
line = trim(line) // "] "
write(line(len_trim(line)+1:), '(I3)') percent
line = trim(line) // "% (" // trim(adjustl(itoa(current))) // "/" // trim(adjustl(itoa(total))) // ")"
if (percent < 100) then
write(error_unit,'(A)') char(13) // char(27) // "[2K" // trim(line) // char(27) // "[1A"
else
write(error_unit,'(A)') char(13) // char(27) // "[2K" // trim(line)
end if
call flush(error_unit)
end subroutine print_progress
pure function itoa(i) result(str)
implicit none
integer(psb_ipk_), intent(in) :: i
character(len=32) :: str
write(str,'(I0)') i
end function itoa
end program main