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 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114
|
!
! Copyright (C) by Argonne National Laboratory
! See COPYRIGHT in top-level directory
!
program main
use mpi
implicit none
integer ierr, i, j, type, count,errs
parameter (count = 4)
integer rank, size, xfersize
integer status(MPI_STATUS_SIZE)
integer blocklens(count), displs(count)
double precision,dimension(:,:),allocatable :: sndbuf, rcvbuf
logical verbose
verbose = .false.
call mtest_init ( ierr )
call mpi_comm_size( MPI_COMM_WORLD, size, ierr )
call mpi_comm_rank( MPI_COMM_WORLD, rank, ierr )
if (size .lt. 2) then
print *, "Must have at least 2 processes"
call MPI_Abort( MPI_COMM_WORLD, 1, ierr )
endif
errs = 0
allocate(sndbuf(7,100))
allocate(rcvbuf(7,100))
do j=1,100
do i=1,7
sndbuf(i,j) = (i+j) * 1.0
enddo
enddo
do i=1,count
blocklens(i) = 7
enddo
! bug occurs when first two displacements are 0
displs(1) = 0
displs(2) = 0
displs(3) = 10
displs(4) = 10
call mpi_type_indexed( count, blocklens, displs*blocklens(1), &
& MPI_DOUBLE_PRECISION, type, ierr )
call mpi_type_commit( type, ierr )
! send using this new type
if (rank .eq. 0) then
call mpi_send( sndbuf(1,1), 1, type, 1, 0, MPI_COMM_WORLD,ierr )
else if (rank .eq. 1) then
xfersize=count * blocklens(1)
call mpi_recv( rcvbuf(1,1), xfersize, MPI_DOUBLE_PRECISION, 0, 0, &
& MPI_COMM_WORLD,status, ierr )
! Values that should be sent
if (verbose) then
! displacement = 0
j=1
do i=1, 7
print*,'sndbuf(',i,j,') = ',sndbuf(i,j)
enddo
! displacement = 10
j=11
do i=1,7
print*,'sndbuf(',i,j,') = ',sndbuf(i,j)
enddo
print*,' '
! Values received
do j=1,count
do i=1,7
print*,'rcvbuf(',i,j,') = ',rcvbuf(i,j)
enddo
enddo
endif
! Error checking
do j=1,2
do i=1,7
if (rcvbuf(i,j) .ne. sndbuf(i,1)) then
print*,'ERROR in rcvbuf(',i,j,')'
print*,'Received ', rcvbuf(i,j),' expected ',sndbuf(i,11)
errs = errs+1
endif
enddo
enddo
do j=3,4
do i=1,7
if (rcvbuf(i,j) .ne. sndbuf(i,11)) then
print*,'ERROR in rcvbuf(',i,j,')'
print*,'Received ', rcvbuf(i,j),' expected ',sndbuf(i,11)
errs = errs+1
endif
enddo
enddo
endif
!
deallocate(sndbuf, rcvbuf)
call mpi_type_free( type, ierr )
call mtest_finalize( errs )
end
|