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 Accumulate Test
9 * Author: James Dinan <dinan@mcs.anl.gov>
10 * Date : December, 2010
12 * This code performs N accumulates into a 2d patch of a shared array. The
13 * array has dimensions [X, Y] and the subarray has dimensions [SUB_X, SUB_Y]
14 * and begins at index [0, 0]. The input and output buffers are specified
15 * using an MPI indexed type.
31 int main(int argc, char **argv) {
32 int itr, i, j, rank, nranks, peer, bufsize, errors;
33 double *win_buf, *src_buf;
36 MTest_Init(&argc, &argv);
38 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
39 MPI_Comm_size(MPI_COMM_WORLD, &nranks);
41 bufsize = XDIM * YDIM * sizeof(double);
42 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &win_buf);
43 MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &src_buf);
45 for (i = 0; i < XDIM*YDIM; i++) {
46 *(win_buf + i) = -1.0;
47 *(src_buf + i) = 1.0 + rank;
50 MPI_Win_create(win_buf, bufsize, 1, MPI_INFO_NULL, MPI_COMM_WORLD, &buf_win);
52 peer = (rank+1) % nranks;
54 /* Perform ITERATIONS strided accumulate operations */
56 for (itr = 0; itr < ITERATIONS; itr++) {
57 MPI_Aint idx_loc[SUB_YDIM];
58 int idx_rem[SUB_YDIM];
59 int blk_len[SUB_YDIM];
60 MPI_Datatype src_type, dst_type;
62 for (i = 0; i < SUB_YDIM; i++) {
63 MPI_Get_address(&src_buf[i*XDIM], &idx_loc[i]);
65 blk_len[i] = SUB_XDIM;
69 MPI_Type_hindexed(SUB_YDIM, blk_len, idx_loc, MPI_DOUBLE, &src_type);
71 MPI_Type_indexed(SUB_YDIM, blk_len, idx_rem, MPI_DOUBLE, &src_type);
73 MPI_Type_indexed(SUB_YDIM, blk_len, idx_rem, MPI_DOUBLE, &dst_type);
75 MPI_Type_commit(&src_type);
76 MPI_Type_commit(&dst_type);
78 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
81 MPI_Accumulate(MPI_BOTTOM, 1, src_type, peer, 0, 1, dst_type, MPI_SUM, buf_win);
83 MPI_Accumulate(src_buf, 1, src_type, peer, 0, 1, dst_type, MPI_SUM, buf_win);
86 MPI_Win_unlock(peer, buf_win);
88 MPI_Type_free(&src_type);
89 MPI_Type_free(&dst_type);
92 MPI_Barrier(MPI_COMM_WORLD);
94 /* Verify that the results are correct */
96 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, rank, 0, buf_win);
98 for (i = 0; i < SUB_XDIM; i++) {
99 for (j = 0; j < SUB_YDIM; j++) {
100 const double actual = *(win_buf + i + j*XDIM);
101 const double expected = -1.0 + (1.0 + ((rank+nranks-1)%nranks)) * (ITERATIONS);
102 if (fabs(actual - expected) > 1.0e-10) {
103 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
104 rank, j, i, expected, actual); );
110 for (i = SUB_XDIM; i < XDIM; i++) {
111 for (j = 0; j < SUB_YDIM; j++) {
112 const double actual = *(win_buf + i + j*XDIM);
113 const double expected = -1.0;
114 if (fabs(actual - expected) > 1.0e-10) {
115 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
116 rank, j, i, expected, actual); );
122 for (i = 0; i < XDIM; i++) {
123 for (j = SUB_YDIM; j < YDIM; j++) {
124 const double actual = *(win_buf + i + j*XDIM);
125 const double expected = -1.0;
126 if (fabs(actual - expected) > 1.0e-10) {
127 SQUELCH( printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
128 rank, j, i, expected, actual); );
134 MPI_Win_unlock(rank, buf_win);
136 MPI_Win_free(&buf_win);
137 MPI_Free_mem(win_buf);
138 MPI_Free_mem(src_buf);
140 MTest_Finalize( errors );
142 return MTestReturnValue( errors );