1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
3 * (C) 2001 by Argonne National Laboratory.
4 * See COPYRIGHT in top-level directory.
7 /* One-Sided MPI 2-D Strided Put Test
9 * Author: James Dinan <dinan@mcs.anl.gov>
12 * This code performs N strided put operations into a 2d patch of a shared
13 * array. The array has dimensions [X, Y] and the subarray has dimensions
14 * [SUB_X, SUB_Y] and begins at index [0, 0]. The input and output buffers are
15 * specified using an MPI datatype.
17 * This test generates a datatype that is relative to MPI_BOTTOM and tests the
18 * RMA implementation's ability to perform the correct transfer.
33 int main(int argc, char **argv) {
34 int i, j, rank, nranks, peer, bufsize, errors;
35 double *win_buf, *src_buf, *dst_buf;
38 MTest_Init(&argc, &argv);
40 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
41 MPI_Comm_size(MPI_COMM_WORLD, &nranks);
43 bufsize = XDIM * YDIM * sizeof(double);
44 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &win_buf);
45 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &src_buf);
46 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &dst_buf);
48 for (i = 0; i < XDIM*YDIM; i++) {
49 *(win_buf + i) = 1.0 + rank;
50 *(src_buf + i) = 1.0 + rank;
53 MPI_Win_create(win_buf, bufsize, 1, MPI_INFO_NULL, MPI_COMM_WORLD, &buf_win);
55 peer = (rank+1) % nranks;
57 /* Perform ITERATIONS strided put operations */
59 for (i = 0; i < ITERATIONS; i++) {
60 MPI_Aint idx_loc[SUB_YDIM];
61 int idx_rem[SUB_YDIM];
62 int blk_len[SUB_YDIM];
63 MPI_Datatype src_type, dst_type;
65 for (j = 0; j < SUB_YDIM; j++) {
66 MPI_Get_address(&src_buf[j*XDIM], &idx_loc[j]);
67 idx_rem[j] = j*XDIM*sizeof(double);
68 blk_len[j] = SUB_XDIM*sizeof(double);
71 MPI_Type_create_hindexed(SUB_YDIM, blk_len, idx_loc, MPI_BYTE, &src_type);
72 MPI_Type_create_indexed_block(SUB_YDIM, SUB_XDIM*sizeof(double), idx_rem, MPI_BYTE, &dst_type);
74 MPI_Type_commit(&src_type);
75 MPI_Type_commit(&dst_type);
77 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
78 MPI_Put(MPI_BOTTOM, 1, src_type, peer, 0, 1, dst_type, buf_win);
79 MPI_Win_unlock(peer, buf_win);
81 MPI_Type_free(&src_type);
82 MPI_Type_free(&dst_type);
85 MPI_Barrier(MPI_COMM_WORLD);
87 /* Verify that the results are correct */
89 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, rank, 0, buf_win);
91 for (i = 0; i < SUB_XDIM; i++) {
92 for (j = 0; j < SUB_YDIM; j++) {
93 const double actual = *(win_buf + i + j*XDIM);
94 const double expected = (1.0 + ((rank+nranks-1)%nranks));
95 if (actual - expected > 1e-10) {
96 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
97 rank, j, i, expected, actual); );
103 for (i = SUB_XDIM; i < XDIM; i++) {
104 for (j = 0; j < SUB_YDIM; j++) {
105 const double actual = *(win_buf + i + j*XDIM);
106 const double expected = 1.0 + rank;
107 if (actual - expected > 1e-10) {
108 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
109 rank, j, i, expected, actual); );
115 for (i = 0; i < XDIM; i++) {
116 for (j = SUB_YDIM; j < YDIM; j++) {
117 const double actual = *(win_buf + i + j*XDIM);
118 const double expected = 1.0 + rank;
119 if (actual - expected > 1e-10) {
120 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
121 rank, j, i, expected, actual); );
127 MPI_Win_unlock(rank, buf_win);
129 MPI_Win_free(&buf_win);
130 MPI_Free_mem(win_buf);
131 MPI_Free_mem(src_buf);
132 MPI_Free_mem(dst_buf);
134 MTest_Finalize( errors );