Skip to content

add test program for MPI_Waitall #97

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merged
merged 1 commit into from
Apr 3, 2025
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion src/mpi.f90
Original file line number Diff line number Diff line change
Expand Up @@ -670,7 +670,7 @@ subroutine MPI_Waitall_proc(count, array_of_requests, array_of_statuses, ierror)
use mpi_c_bindings, only: c_mpi_waitall
integer, intent(in) :: count
integer, intent(inout) :: array_of_requests(count)
integer :: array_of_statuses(*)
integer, intent(out) :: array_of_statuses(*)
integer, optional, intent(out) :: ierror
call c_mpi_waitall(count, array_of_requests, array_of_statuses, ierror)
end subroutine
Expand Down
67 changes: 67 additions & 0 deletions tests/waitall_1.f90
Original file line number Diff line number Diff line change
@@ -0,0 +1,67 @@
program waitall_1
use mpi
implicit none

integer, parameter :: lbuf = 100
integer, parameter :: ntype_real = MPI_DOUBLE_PRECISION
real(8), dimension(lbuf,2) :: sbuf21, sbuf22
real(8), dimension(lbuf,2) :: rbuf21, rbuf22
integer, dimension(4) :: reqs
integer :: ierr, myrank, nprocs, iproc_tp, iproc_tm
integer :: comm_all, tag
real(8) :: expected_rbuf21, expected_rbuf22
real(8), parameter :: tol = 1.0e-10

call MPI_Init(ierr)
comm_all = MPI_COMM_WORLD
call MPI_Comm_rank(comm_all, myrank, ierr)
call MPI_Comm_size(comm_all, nprocs, ierr)

tag = 1

iproc_tp = mod(myrank + 1, nprocs)
iproc_tm = mod(myrank - 1 + nprocs, nprocs)

sbuf21 = myrank + 0.1
sbuf22 = myrank + 0.2

call MPI_Isend(sbuf21, lbuf*2, ntype_real, iproc_tp, tag, &
comm_all, reqs(1), ierr)

call MPI_Isend(sbuf22, lbuf*2, ntype_real, iproc_tm, tag, &
comm_all, reqs(2), ierr)

call MPI_Irecv(rbuf21, lbuf*2, ntype_real, iproc_tm, tag, &
comm_all, reqs(3), ierr)

call MPI_Irecv(rbuf22, lbuf*2, ntype_real, iproc_tp, tag, &
comm_all, reqs(4), ierr)

call MPI_Waitall(4, reqs, MPI_STATUSES_IGNORE, ierr)

expected_rbuf21 = real(iproc_tm) + 0.1
expected_rbuf22 = real(iproc_tp) + 0.2

if (any(abs(rbuf21 - expected_rbuf21) > tol)) then
write(*,*) 'Rank', myrank, ': rbuf21 validation failed'
write(*,*) 'Expected:', expected_rbuf21
write(*,*) 'Received:', rbuf21(1,1)
error stop 'MPI communication error in rbuf21'
end if

if (any(abs(rbuf22 - expected_rbuf22) > tol)) then
write(*,*) 'Rank', myrank, ': rbuf22 validation failed'
write(*,*) 'Expected:', expected_rbuf22
write(*,*) 'Received:', rbuf22(1,1)
error stop 'MPI communication error in rbuf22'
end if

if (myrank == 0) then
write(*,*) 'Process', myrank, 'received data:'
write(*,*) 'rbuf21(1,1) =', rbuf21(1,1), ' (Expected:', expected_rbuf21, ')'
write(*,*) 'rbuf22(1,1) =', rbuf22(1,1), ' (Expected:', expected_rbuf22, ')'
end if

call MPI_Finalize(ierr)

end program waitall_1