1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88
|
! This file created from f77/pt2pt/sendrecvreplf.f with f77tof90
!
! Copyright (C) by Argonne National Laboratory
! See COPYRIGHT in top-level directory
!
! This program is based on the allpair.f test from the MPICH-1 test
! (test/pt2pt/allpair.f), which in turn was inspired by a bug report from
! fsset@corelli.lerc.nasa.gov (Scott Townsend)
program sendrecvrepl
use mpi
integer ierr, errs, comm
logical mtestGetIntraComm
logical verbose
common /flags/ verbose
errs = 0
verbose = .false.
! verbose = .true.
call MTest_Init( ierr )
do while ( mtestGetIntraComm( comm, 2, .false. ) )
call test_pair_sendrecvrepl( comm, errs )
call mtestFreeComm( comm )
enddo
!
call MTest_Finalize( errs )
!
end
!
subroutine test_pair_sendrecvrepl( comm, errs )
use mpi
integer comm, errs
integer rank, size, ierr, next, prev, tag, count, i
integer TEST_SIZE
parameter (TEST_SIZE=2000)
integer status(MPI_STATUS_SIZE)
real send_buf(TEST_SIZE), recv_buf(TEST_SIZE)
logical verbose
common /flags/ verbose
!
if (verbose) then
print *, ' Sendrecv replace'
endif
!
call mpi_comm_rank( comm, rank, ierr )
call mpi_comm_size( comm, size, ierr )
next = rank + 1
if (next .ge. size) next = 0
!
prev = rank - 1
if (prev .lt. 0) prev = size - 1
!
tag = 4456
count = TEST_SIZE / 3
if (rank .eq. 0) then
!
call init_test_data(recv_buf, TEST_SIZE)
!
do 11 i = count+1,TEST_SIZE
recv_buf(i) = 0.0
11 continue
!
call MPI_Sendrecv_replace(recv_buf, count, MPI_REAL, &
& next, tag, next, tag, &
& comm, status, ierr)
call msg_check( recv_buf, next, tag, count, status, TEST_SIZE, &
& 'sendrecvreplace', errs )
else if (prev .eq. 0) then
call clear_test_data(recv_buf,TEST_SIZE)
call MPI_Recv(recv_buf, TEST_SIZE, MPI_REAL, &
& MPI_ANY_SOURCE, MPI_ANY_TAG, comm, &
& status, ierr)
call msg_check( recv_buf, prev, tag, count, status, TEST_SIZE, &
& 'recv/send for replace', errs )
call MPI_Send(recv_buf, count, MPI_REAL, prev, tag, &
& comm, ierr)
end if
!
end
|