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 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176
|
C
C Copyright (C) by Argonne National Laboratory
C See COPYRIGHT in top-level directory
C
program main
implicit none
include 'mpif.h'
integer errs, ierr, i, intsize
integer type1, type2, type3, type4, type5
integer max_asizev
parameter (max_asizev = 10)
include 'typeaints.h'
integer blocklens(max_asizev), dtypes(max_asizev)
integer displs(max_asizev)
integer recvbuf(6*max_asizev)
integer sendbuf(max_asizev), status(MPI_STATUS_SIZE)
integer rank, size
errs = 0
call mtest_init( ierr )
call mpi_comm_size( MPI_COMM_WORLD, size, ierr )
call mpi_comm_rank( MPI_COMM_WORLD, rank, ierr )
C
call mpi_type_size( MPI_INTEGER, intsize, ierr )
C
aintv(1) = 0
aintv(2) = 3 * intsize
call mpi_type_create_resized( MPI_INTEGER, aintv(1), aintv(2),
& type1, ierr )
call mpi_type_commit( type1, ierr )
aintv(1) = -1
aintv(2) = -1
call mpi_type_get_extent( type1, aintv(1), aintv(2), ierr )
if (aintv(1) .ne. 0) then
errs = errs + 1
print *, 'Did not get expected lb'
endif
if (aintv(2) .ne. 3*intsize) then
errs = errs + 1
print *, 'Did not get expected extent'
endif
aintv(1) = -1
aintv(2) = -1
call mpi_type_get_true_extent( type1, aintv(1), aintv(2), ierr )
if (aintv(1) .ne. 0) then
errs = errs + 1
print *, 'Did not get expected true lb'
endif
if (aintv(2) .ne. intsize) then
errs = errs + 1
print *, 'Did not get expected true extent (', aintv(2), ') ',
& ' expected ', intsize
endif
C
do i=1,10
blocklens(i) = 1
aintv(i) = (i-1) * 3 * intsize
enddo
call mpi_type_create_hindexed( 10, blocklens, aintv,
& MPI_INTEGER, type2, ierr )
call mpi_type_commit( type2, ierr )
C
aint = 3 * intsize
call mpi_type_create_hvector( 10, 1, aint, MPI_INTEGER, type3,
& ierr )
call mpi_type_commit( type3, ierr )
C
do i=1,10
blocklens(i) = 1
dtypes(i) = MPI_INTEGER
aintv(i) = (i-1) * 3 * intsize
enddo
call mpi_type_create_struct( 10, blocklens, aintv, dtypes,
& type4, ierr )
call mpi_type_commit( type4, ierr )
do i=1,10
displs(i) = (i-1) * 3
enddo
call mpi_type_create_indexed_block( 10, 1, displs,
& MPI_INTEGER, type5, ierr )
call mpi_type_commit( type5, ierr )
C
C Using each time, send and receive using these types
do i=1, max_asizev*3
recvbuf(i) = -1
enddo
do i=1, max_asizev
sendbuf(i) = i
enddo
call mpi_sendrecv( sendbuf, max_asizev, MPI_INTEGER, rank, 0,
& recvbuf, max_asizev, type1, rank, 0,
& MPI_COMM_WORLD, status, ierr )
do i=1, max_asizev
if (recvbuf(1+(i-1)*3) .ne. i ) then
errs = errs + 1
print *, 'type1:', i, 'th element = ', recvbuf(1+(i-1)*3)
endif
enddo
C
do i=1, max_asizev*3
recvbuf(i) = -1
enddo
do i=1, max_asizev
sendbuf(i) = i
enddo
call mpi_sendrecv( sendbuf, max_asizev, MPI_INTEGER, rank, 0,
& recvbuf, 1, type2, rank, 0,
& MPI_COMM_WORLD, status, ierr )
do i=1, max_asizev
if (recvbuf(1+(i-1)*3) .ne. i ) then
errs = errs + 1
print *, 'type2:', i, 'th element = ', recvbuf(1+(i-1)*3)
endif
enddo
C
do i=1, max_asizev*3
recvbuf(i) = -1
enddo
do i=1, max_asizev
sendbuf(i) = i
enddo
call mpi_sendrecv( sendbuf, max_asizev, MPI_INTEGER, rank, 0,
& recvbuf, 1, type3, rank, 0,
& MPI_COMM_WORLD, status, ierr )
do i=1, max_asizev
if (recvbuf(1+(i-1)*3) .ne. i ) then
errs = errs + 1
print *, 'type3:', i, 'th element = ', recvbuf(1+(i-1)*3)
endif
enddo
C
do i=1, max_asizev*3
recvbuf(i) = -1
enddo
do i=1, max_asizev
sendbuf(i) = i
enddo
call mpi_sendrecv( sendbuf, max_asizev, MPI_INTEGER, rank, 0,
& recvbuf, 1, type4, rank, 0,
& MPI_COMM_WORLD, status, ierr )
do i=1, max_asizev
if (recvbuf(1+(i-1)*3) .ne. i ) then
errs = errs + 1
print *, 'type4:', i, 'th element = ', recvbuf(1+(i-1)*3)
endif
enddo
C
do i=1, max_asizev*3
recvbuf(i) = -1
enddo
do i=1, max_asizev
sendbuf(i) = i
enddo
call mpi_sendrecv( sendbuf, max_asizev, MPI_INTEGER, rank, 0,
& recvbuf, 1, type5, rank, 0,
& MPI_COMM_WORLD, status, ierr )
do i=1, max_asizev
if (recvbuf(1+(i-1)*3) .ne. i ) then
errs = errs + 1
print *, 'type5:', i, 'th element = ', recvbuf(1+(i-1)*3)
endif
enddo
C
call mpi_type_free( type1, ierr )
call mpi_type_free( type2, ierr )
call mpi_type_free( type3, ierr )
call mpi_type_free( type4, ierr )
call mpi_type_free( type5, ierr )
call mtest_finalize( errs )
end
|