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 subarray type.
31 int main(int argc, char **argv) {
32 int 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 (i = 0; i < ITERATIONS; i++) {
58 int src_arr_sizes[2] = { XDIM, YDIM };
59 int src_arr_subsizes[2] = { SUB_XDIM, SUB_YDIM };
60 int src_arr_starts[2] = { 0, 0 };
61 int dst_arr_sizes[2] = { XDIM, YDIM };
62 int dst_arr_subsizes[2] = { SUB_XDIM, SUB_YDIM };
63 int dst_arr_starts[2] = { 0, 0 };
64 MPI_Datatype src_type, dst_type;
66 MPI_Type_create_subarray(ndims, src_arr_sizes, src_arr_subsizes, src_arr_starts,
67 MPI_ORDER_C, MPI_DOUBLE, &src_type);
69 MPI_Type_create_subarray(ndims, dst_arr_sizes, dst_arr_subsizes, dst_arr_starts,
70 MPI_ORDER_C, MPI_DOUBLE, &dst_type);
72 MPI_Type_commit(&src_type);
73 MPI_Type_commit(&dst_type);
75 MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
77 MPI_Accumulate(src_buf, 1, src_type, peer, 0, 1, dst_type, MPI_SUM, 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 + (1.0 + ((rank+nranks-1)%nranks)) * (ITERATIONS);
95 if (fabs(actual - expected) > 1.0e-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;
107 if (fabs(actual - expected) > 1.0e-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;
119 if (fabs(actual - expected) > 1.0e-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);
133 MTest_Finalize( errors );
135 return MTestReturnValue( errors );