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)
35 int i, j, rank, nranks, peer, bufsize, errors;
36 double *win_buf, *src_buf, *dst_buf;
39 MTest_Init(&argc, &argv);
41 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
42 MPI_Comm_size(MPI_COMM_WORLD, &nranks);
44 bufsize = XDIM * YDIM * sizeof(double);
45 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &win_buf);
46 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &src_buf);
47 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &dst_buf);
49 for (i = 0; i < XDIM * YDIM; i++) {
50 *(win_buf + i) = 1.0 + rank;
51 *(src_buf + i) = 1.0 + rank;
54 MPI_Win_create(win_buf, bufsize, 1, MPI_INFO_NULL, MPI_COMM_WORLD, &buf_win);
56 peer = (rank + 1) % nranks;
58 /* Perform ITERATIONS strided put operations */
60 for (i = 0; i < ITERATIONS; i++) {
61 MPI_Aint idx_loc[SUB_YDIM];
62 int idx_rem[SUB_YDIM];
63 int blk_len[SUB_YDIM];
64 MPI_Datatype src_type, dst_type;
66 for (j = 0; j < SUB_YDIM; j++) {
67 MPI_Get_address(&src_buf[j * XDIM], &idx_loc[j]);
68 idx_rem[j] = j * XDIM * sizeof(double);
69 blk_len[j] = SUB_XDIM * sizeof(double);
72 MPI_Type_create_hindexed(SUB_YDIM, blk_len, idx_loc, MPI_BYTE, &src_type);
73 MPI_Type_create_indexed_block(SUB_YDIM, SUB_XDIM * sizeof(double), idx_rem, MPI_BYTE,
76 MPI_Type_commit(&src_type);
77 MPI_Type_commit(&dst_type);
79 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
80 MPI_Put(MPI_BOTTOM, 1, src_type, peer, 0, 1, dst_type, buf_win);
81 MPI_Win_unlock(peer, buf_win);
83 MPI_Type_free(&src_type);
84 MPI_Type_free(&dst_type);
87 MPI_Barrier(MPI_COMM_WORLD);
89 /* Verify that the results are correct */
91 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, rank, 0, buf_win);
93 for (i = 0; i < SUB_XDIM; i++) {
94 for (j = 0; j < SUB_YDIM; j++) {
95 const double actual = *(win_buf + i + j * XDIM);
96 const double expected = (1.0 + ((rank + nranks - 1) % nranks));
97 if (actual - expected > 1e-10) {
98 SQUELCH(printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
99 rank, j, i, expected, actual););
105 for (i = SUB_XDIM; i < XDIM; i++) {
106 for (j = 0; j < SUB_YDIM; j++) {
107 const double actual = *(win_buf + i + j * XDIM);
108 const double expected = 1.0 + rank;
109 if (actual - expected > 1e-10) {
110 SQUELCH(printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
111 rank, j, i, expected, actual););
117 for (i = 0; i < XDIM; i++) {
118 for (j = SUB_YDIM; j < YDIM; j++) {
119 const double actual = *(win_buf + i + j * XDIM);
120 const double expected = 1.0 + rank;
121 if (actual - expected > 1e-10) {
122 SQUELCH(printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
123 rank, j, i, expected, actual););
129 MPI_Win_unlock(rank, buf_win);
131 MPI_Win_free(&buf_win);
132 MPI_Free_mem(win_buf);
133 MPI_Free_mem(src_buf);
134 MPI_Free_mem(dst_buf);
136 MTest_Finalize(errors);