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 an arbitrary base address
18 * in memory and tests the RMA implementation's ability to perform the correct
35 int main(int argc, char **argv)
37 int i, j, rank, nranks, peer, bufsize, errors;
38 double *win_buf, *src_buf, *dst_buf;
41 MTest_Init(&argc, &argv);
43 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
44 MPI_Comm_size(MPI_COMM_WORLD, &nranks);
46 bufsize = XDIM * YDIM * sizeof(double);
47 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &win_buf);
48 /* Alloc_mem is not required for the origin buffers for RMA operations -
49 * just for the Win_create memory */
50 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &src_buf);
51 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &dst_buf);
53 for (i = 0; i < XDIM * YDIM; i++) {
54 *(win_buf + i) = 1.0 + rank;
55 *(src_buf + i) = 1.0 + rank;
58 MPI_Win_create(win_buf, bufsize, 1, MPI_INFO_NULL, MPI_COMM_WORLD, &buf_win);
60 peer = (rank + 1) % nranks;
62 /* Perform ITERATIONS strided put operations */
64 for (i = 0; i < ITERATIONS; i++) {
65 MPI_Aint idx_loc[SUB_YDIM];
66 int idx_rem[SUB_YDIM];
67 int blk_len[SUB_YDIM];
68 MPI_Datatype src_type, dst_type;
70 void *base_ptr = dst_buf;
73 MPI_Get_address(base_ptr, &base_int);
75 for (j = 0; j < SUB_YDIM; j++) {
76 MPI_Get_address(&src_buf[j * XDIM], &idx_loc[j]);
77 idx_loc[j] = idx_loc[j] - base_int;
78 idx_rem[j] = j * XDIM * sizeof(double);
79 blk_len[j] = SUB_XDIM * sizeof(double);
82 MPI_Type_create_hindexed(SUB_YDIM, blk_len, idx_loc, MPI_BYTE, &src_type);
83 MPI_Type_create_indexed_block(SUB_YDIM, SUB_XDIM * sizeof(double), idx_rem, MPI_BYTE,
86 MPI_Type_commit(&src_type);
87 MPI_Type_commit(&dst_type);
89 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
90 MPI_Put(base_ptr, 1, src_type, peer, 0, 1, dst_type, buf_win);
91 MPI_Win_unlock(peer, buf_win);
93 MPI_Type_free(&src_type);
94 MPI_Type_free(&dst_type);
97 MPI_Barrier(MPI_COMM_WORLD);
99 /* Verify that the results are correct */
101 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, rank, 0, buf_win);
103 for (i = 0; i < SUB_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 + nranks - 1) % nranks));
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 = SUB_XDIM; i < XDIM; i++) {
116 for (j = 0; j < SUB_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 for (i = 0; i < XDIM; i++) {
128 for (j = SUB_YDIM; j < YDIM; j++) {
129 const double actual = *(win_buf + i + j * XDIM);
130 const double expected = 1.0 + rank;
131 if (actual - expected > 1e-10) {
132 SQUELCH(printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
133 rank, j, i, expected, actual););
139 MPI_Win_unlock(rank, buf_win);
141 MPI_Win_free(&buf_win);
142 MPI_Free_mem(win_buf);
143 MPI_Free_mem(src_buf);
144 MPI_Free_mem(dst_buf);
146 MTest_Finalize(errors);