1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
3 * (C) 2012 by Argonne National Laboratory.
4 * See COPYRIGHT in top-level directory.
7 /* One-Sided MPI 2-D Strided Accumulate Test
9 * Author: James Dinan <dinan@mcs.anl.gov>
10 * Date : November, 2012
12 * This code performs N strided put operations followed by get operations into
13 * a 2d patch of a shared array. The array has dimensions [X, Y] and the
14 * subarray has dimensions [SUB_X, SUB_Y] and begins at index [0, 0]. The
15 * input and output buffers are specified using an MPI indexed type.
31 int main(int argc, char **argv) {
32 int rank, nranks, rank_world, nranks_world;
33 int i, j, peer, bufsize, errors;
34 double *win_buf, *src_buf, *dst_buf;
38 MTest_Init(&argc, &argv);
40 MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
41 MPI_Comm_size(MPI_COMM_WORLD, &nranks_world);
43 MPI_Comm_split_type(MPI_COMM_WORLD, MPI_COMM_TYPE_SHARED, rank, MPI_INFO_NULL, &shr_comm);
45 MPI_Comm_rank(shr_comm, &rank);
46 MPI_Comm_size(shr_comm, &nranks);
48 bufsize = XDIM * YDIM * sizeof(double);
49 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &src_buf);
50 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &dst_buf);
52 MPI_Win_allocate_shared(bufsize, 1, MPI_INFO_NULL, shr_comm, &win_buf, &buf_win);
54 MPI_Win_fence(0, buf_win);
56 for (i = 0; i < XDIM*YDIM; i++) {
57 *(win_buf + i) = -1.0;
58 *(src_buf + i) = 1.0 + rank;
61 MPI_Win_fence(0, buf_win);
63 peer = (rank+1) % nranks;
65 /* Perform ITERATIONS strided accumulate operations */
67 for (i = 0; i < ITERATIONS; i++) {
68 int idx_rem[SUB_YDIM];
69 int blk_len[SUB_YDIM];
70 MPI_Datatype src_type, dst_type;
72 for (j = 0; j < SUB_YDIM; j++) {
74 blk_len[j] = SUB_XDIM;
77 MPI_Type_indexed(SUB_YDIM, blk_len, idx_rem, MPI_DOUBLE, &src_type);
78 MPI_Type_indexed(SUB_YDIM, blk_len, idx_rem, MPI_DOUBLE, &dst_type);
80 MPI_Type_commit(&src_type);
81 MPI_Type_commit(&dst_type);
83 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
84 MPI_Put(src_buf, 1, src_type, peer, 0, 1, dst_type, buf_win);
85 MPI_Win_unlock(peer, buf_win);
87 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
88 MPI_Get(dst_buf, 1, src_type, peer, 0, 1, dst_type, buf_win);
89 MPI_Win_unlock(peer, buf_win);
91 MPI_Type_free(&src_type);
92 MPI_Type_free(&dst_type);
95 MPI_Barrier(shr_comm);
97 /* Verify that the results are correct */
99 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, rank, 0, buf_win);
101 for (i = 0; i < SUB_XDIM; i++) {
102 for (j = 0; j < SUB_YDIM; j++) {
103 const double actual = *(win_buf + i + j*XDIM);
104 const double expected = (1.0 + ((rank+nranks-1)%nranks));
105 if (fabs(actual - expected) > 1.0e-10) {
106 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
107 rank, j, i, expected, actual); );
113 for (i = SUB_XDIM; i < XDIM; i++) {
114 for (j = 0; j < SUB_YDIM; j++) {
115 const double actual = *(win_buf + i + j*XDIM);
116 const double expected = -1.0;
117 if (fabs(actual - expected) > 1.0e-10) {
118 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
119 rank, j, i, expected, actual); );
125 for (i = 0; i < XDIM; i++) {
126 for (j = SUB_YDIM; j < YDIM; j++) {
127 const double actual = *(win_buf + i + j*XDIM);
128 const double expected = -1.0;
129 if (fabs(actual - expected) > 1.0e-10) {
130 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
131 rank, j, i, expected, actual); );
137 MPI_Win_unlock(rank, buf_win);
139 MPI_Win_free(&buf_win);
140 MPI_Free_mem(src_buf);
141 MPI_Free_mem(dst_buf);
142 MPI_Comm_free(&shr_comm);
144 MTest_Finalize( errors );
146 return MTestReturnValue( errors );