Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
fix+activate rma test
[simgrid.git] / teshsuite / smpi / mpich3-test / rma / strided_acc_onelock.c
1 /* -*- Mode: C; c-basic-offset:4 ; indent-tabs-mode:nil ; -*- */
2 /*
3  *  (C) 2001 by Argonne National Laboratory.
4  *      See COPYRIGHT in top-level directory.
5  */
6
7 /* One-Sided MPI 2-D Strided Accumulate Test
8  *
9  * Author: James Dinan <dinan@mcs.anl.gov>
10  * Date  : December, 2010
11  *
12  * This code performs one-sided accumulate into a 2d patch of a shared array.
13  */
14
15 #include <stdio.h>
16 #include <stdlib.h>
17 #include <math.h>
18 #include <mpi.h>
19 #include "mpitest.h"
20 #include "squelch.h"
21
22 #define XDIM 1024
23 #define YDIM 1024
24 #define ITERATIONS 3
25
26 int main(int argc, char **argv)
27 {
28     int i, j, rank, nranks, peer, bufsize, errors;
29     double *buffer, *src_buf;
30     MPI_Win buf_win;
31
32     MTest_Init(&argc, &argv);
33
34     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
35     MPI_Comm_size(MPI_COMM_WORLD, &nranks);
36
37     bufsize = XDIM * YDIM * sizeof(double);
38     MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &buffer);
39     MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &src_buf);
40
41     for (i = 0; i < XDIM * YDIM; i++) {
42         *(buffer + i) = 1.0 + rank;
43         *(src_buf + i) = 1.0 + rank;
44     }
45
46     MPI_Win_create(buffer, bufsize, 1, MPI_INFO_NULL, MPI_COMM_WORLD, &buf_win);
47
48     peer = (rank + 1) % nranks;
49
50     for (i = 0; i < ITERATIONS; i++) {
51
52         MPI_Win_lock(MPI_LOCK_EXCLUSIVE, peer, 0, buf_win);
53
54         for (j = 0; j < YDIM; j++) {
55             MPI_Accumulate(src_buf + j * XDIM, XDIM, MPI_DOUBLE, peer,
56                            j * XDIM * sizeof(double), XDIM, MPI_DOUBLE, MPI_SUM, buf_win);
57         }
58
59         MPI_Win_unlock(peer, buf_win);
60     }
61
62     MPI_Barrier(MPI_COMM_WORLD);
63
64     MPI_Win_lock(MPI_LOCK_EXCLUSIVE, rank, 0, buf_win);
65     for (i = errors = 0; i < XDIM; i++) {
66         for (j = 0; j < YDIM; j++) {
67             const double actual = *(buffer + i + j * XDIM);
68             const double expected =
69                 (1.0 + rank) + (1.0 + ((rank + nranks - 1) % nranks)) * (ITERATIONS);
70             if (fabs(actual - expected) > 1.0e-10) {
71                 SQUELCH(printf("%d: Data validation failed at [%d, %d] expected=%f actual=%f\n",
72                                rank, j, i, expected, actual););
73                 errors++;
74                 fflush(stdout);
75             }
76         }
77     }
78     MPI_Win_unlock(rank, buf_win);
79
80     MPI_Win_free(&buf_win);
81     MPI_Free_mem(buffer);
82     MPI_Free_mem(src_buf);
83
84     MTest_Finalize(errors);
85     MPI_Finalize();
86     return MTestReturnValue(errors);
87 }