1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150
|
/*
* Copyright (C) by Argonne National Laboratory
* See COPYRIGHT in top-level directory
*/
/* One-Sided MPI 2-D Strided Accumulate Test
*
* Author: James Dinan <dinan@mcs.anl.gov>
* Date : November, 2012
*
* This code performs N strided put operations followed by get operations into
* a 2d patch of a shared array. The array has dimensions [X, Y] and the
* subarray has dimensions [SUB_X, SUB_Y] and begins at index [0, 0]. The
* input and output buffers are specified using an MPI indexed type.
*/
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <mpi.h>
#include "mpitest.h"
#include "squelch.h"
#define XDIM 8
#define YDIM 1024
#define SUB_XDIM 1
#define SUB_YDIM 2
#define ITERATIONS 10
int main(int argc, char **argv)
{
int rank, nranks, rank_world, nranks_world;
int i, j, peer, bufsize, errors;
double *win_buf, *src_buf, *dst_buf;
MPI_Win buf_win;
MPI_Comm shr_comm;
MTest_Init(&argc, &argv);
MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
MPI_Comm_size(MPI_COMM_WORLD, &nranks_world);
MPI_Comm_split_type(MPI_COMM_WORLD, MPI_COMM_TYPE_SHARED, rank_world, MPI_INFO_NULL, &shr_comm);
MPI_Comm_rank(shr_comm, &rank);
MPI_Comm_size(shr_comm, &nranks);
bufsize = XDIM * YDIM * sizeof(double);
MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &src_buf);
MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &dst_buf);
MPI_Win_allocate_shared(bufsize, 1, MPI_INFO_NULL, shr_comm, &win_buf, &buf_win);
MPI_Win_fence(0, buf_win);
for (i = 0; i < XDIM * YDIM; i++) {
*(win_buf + i) = -1.0;
*(src_buf + i) = 1.0 + rank;
}
MPI_Win_fence(0, buf_win);
peer = (rank + 1) % nranks;
/* Perform ITERATIONS strided accumulate operations */
for (i = 0; i < ITERATIONS; i++) {
int idx_rem[SUB_YDIM];
int blk_len[SUB_YDIM];
MPI_Datatype src_type, dst_type;
for (j = 0; j < SUB_YDIM; j++) {
idx_rem[j] = j * XDIM;
blk_len[j] = SUB_XDIM;
}
MPI_Type_indexed(SUB_YDIM, blk_len, idx_rem, MPI_DOUBLE, &src_type);
MPI_Type_indexed(SUB_YDIM, blk_len, idx_rem, MPI_DOUBLE, &dst_type);
MPI_Type_commit(&src_type);
MPI_Type_commit(&dst_type);
/* PUT */
MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
MPI_Get_accumulate(src_buf, 1, src_type, dst_buf, 1, src_type, peer, 0,
1, dst_type, MPI_REPLACE, buf_win);
MPI_Win_unlock(peer, buf_win);
/* GET */
MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
MPI_Get_accumulate(src_buf, 1, src_type, dst_buf, 1, src_type, peer, 0,
1, dst_type, MPI_NO_OP, buf_win);
MPI_Win_unlock(peer, buf_win);
MPI_Type_free(&src_type);
MPI_Type_free(&dst_type);
}
MPI_Barrier(MPI_COMM_WORLD);
/* Verify that the results are correct */
MPI_Win_lock(MPI_LOCK_EXCLUSIVE, rank, 0, buf_win);
errors = 0;
for (i = 0; i < SUB_XDIM; i++) {
for (j = 0; j < SUB_YDIM; j++) {
const double actual = *(win_buf + i + j * XDIM);
const double expected = (1.0 + ((rank + nranks - 1) % nranks));
if (fabs(actual - expected) > 1.0e-10) {
SQUELCH(printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
rank, j, i, expected, actual););
errors++;
fflush(stdout);
}
}
}
for (i = SUB_XDIM; i < XDIM; i++) {
for (j = 0; j < SUB_YDIM; j++) {
const double actual = *(win_buf + i + j * XDIM);
const double expected = -1.0;
if (fabs(actual - expected) > 1.0e-10) {
SQUELCH(printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
rank, j, i, expected, actual););
errors++;
fflush(stdout);
}
}
}
for (i = 0; i < XDIM; i++) {
for (j = SUB_YDIM; j < YDIM; j++) {
const double actual = *(win_buf + i + j * XDIM);
const double expected = -1.0;
if (fabs(actual - expected) > 1.0e-10) {
SQUELCH(printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
rank, j, i, expected, actual););
errors++;
fflush(stdout);
}
}
}
MPI_Win_unlock(rank, buf_win);
MPI_Win_free(&buf_win);
MPI_Free_mem(src_buf);
MPI_Free_mem(dst_buf);
MPI_Comm_free(&shr_comm);
MTest_Finalize(errors);
return MTestReturnValue(errors);
}
|