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
|
/*
* Copyright (C) by Argonne National Laboratory
* See COPYRIGHT in top-level directory
*/
#include <stdio.h>
#include <mpi.h>
#include "mpitest.h"
#define PUT_SIZE 1
#define GET_SIZE 100000
#define WIN_SIZE (PUT_SIZE+GET_SIZE)
#define LOOP 100
int win_buf[WIN_SIZE], orig_get_buf[GET_SIZE], orig_put_buf[PUT_SIZE];
int main(int argc, char **argv)
{
MPI_Win win;
int i, k, rank, nproc;
int orig_rank = 0, dest_rank = 1;
int errs = 0;
MPI_Group comm_group, orig_group, dest_group;
MTest_Init(&argc, &argv);
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &nproc);
MPI_Comm_group(MPI_COMM_WORLD, &comm_group);
MPI_Group_incl(comm_group, 1, &orig_rank, &orig_group);
MPI_Group_incl(comm_group, 1, &dest_rank, &dest_group);
for (i = 0; i < PUT_SIZE; i++)
orig_put_buf[i] = 1;
for (i = 0; i < GET_SIZE; i++)
orig_get_buf[i] = 0;
for (i = 0; i < WIN_SIZE; i++)
win_buf[i] = 1;
MPI_Win_create(win_buf, sizeof(int) * WIN_SIZE, sizeof(int), MPI_INFO_NULL,
MPI_COMM_WORLD, &win);
for (k = 0; k < LOOP; k++) {
/* test for FENCE */
if (rank == orig_rank) {
MPI_Win_fence(MPI_MODE_NOPRECEDE, win);
MPI_Get(orig_get_buf, GET_SIZE, MPI_INT,
dest_rank, PUT_SIZE /*disp */ , GET_SIZE, MPI_INT, win);
MPI_Put(orig_put_buf, PUT_SIZE, MPI_INT,
dest_rank, 0 /*disp */ , PUT_SIZE, MPI_INT, win);
MPI_Win_fence(MPI_MODE_NOSUCCEED, win);
/* check GET result values */
for (i = 0; i < GET_SIZE; i++) {
if (orig_get_buf[i] != 1) {
printf("LOOP=%d, FENCE, orig_get_buf[%d] = %d, expected 1.\n",
k, i, orig_get_buf[i]);
errs++;
}
}
} else if (rank == dest_rank) {
MPI_Win_fence(MPI_MODE_NOPRECEDE, win);
MPI_Win_fence(MPI_MODE_NOSUCCEED, win);
/* modify the last element in window */
MPI_Win_lock(MPI_LOCK_SHARED, rank, 0, win);
win_buf[WIN_SIZE - 1] = 2;
MPI_Win_unlock(rank, win);
}
MPI_Barrier(MPI_COMM_WORLD);
/* reset buffers */
for (i = 0; i < PUT_SIZE; i++)
orig_put_buf[i] = 1;
for (i = 0; i < GET_SIZE; i++)
orig_get_buf[i] = 0;
MPI_Win_lock(MPI_LOCK_SHARED, rank, 0, win);
for (i = 0; i < WIN_SIZE; i++)
win_buf[i] = 1;
MPI_Win_unlock(rank, win);
MPI_Barrier(MPI_COMM_WORLD);
/* test for PSCW */
if (rank == orig_rank) {
MPI_Win_start(dest_group, 0, win);
MPI_Get(orig_get_buf, GET_SIZE, MPI_INT,
dest_rank, PUT_SIZE /*disp */ , GET_SIZE, MPI_INT, win);
MPI_Put(orig_put_buf, PUT_SIZE, MPI_INT,
dest_rank, 0 /*disp */ , PUT_SIZE, MPI_INT, win);
MPI_Win_complete(win);
/* check GET result values */
for (i = 0; i < GET_SIZE; i++) {
if (orig_get_buf[i] != 1) {
printf("LOOP=%d, PSCW, orig_get_buf[%d] = %d, expected 1.\n",
k, i, orig_get_buf[i]);
errs++;
}
}
} else if (rank == dest_rank) {
MPI_Win_post(orig_group, 0, win);
MPI_Win_wait(win);
/* modify the last element in window */
MPI_Win_lock(MPI_LOCK_SHARED, rank, 0, win);
win_buf[WIN_SIZE - 1] = 2;
MPI_Win_unlock(rank, win);
}
MPI_Barrier(MPI_COMM_WORLD);
/* reset buffers */
for (i = 0; i < PUT_SIZE; i++)
orig_put_buf[i] = 1;
for (i = 0; i < GET_SIZE; i++)
orig_get_buf[i] = 0;
MPI_Win_lock(MPI_LOCK_SHARED, rank, 0, win);
for (i = 0; i < WIN_SIZE; i++)
win_buf[i] = 1;
MPI_Win_unlock(rank, win);
MPI_Barrier(MPI_COMM_WORLD);
}
MPI_Win_free(&win);
MPI_Group_free(&orig_group);
MPI_Group_free(&dest_group);
MPI_Group_free(&comm_group);
MTest_Finalize(errs);
return MTestReturnValue(errs);
}
|