File: scatterv.c

package info (click to toggle)
mpich 4.3.0%2Breally4.2.1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 419,120 kB
  • sloc: ansic: 1,215,557; cpp: 74,755; javascript: 40,763; f90: 20,649; sh: 18,463; xml: 14,418; python: 14,397; perl: 13,772; makefile: 9,279; fortran: 8,063; java: 4,553; asm: 324; ruby: 176; lisp: 19; php: 8; sed: 4
file content (175 lines) | stat: -rw-r--r-- 5,487 bytes parent folder | download | duplicates (4)
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
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
/*
 * Copyright (C) by Argonne National Laboratory
 *     See COPYRIGHT in top-level directory
 */

#include "mpi.h"
#include <stdlib.h>
#include <stdio.h>
#include "mpitest.h"

/* Prototypes for picky compilers */
void SetData(double *, double *, int, int, int, int, int, int);
int CheckData(double *, int, int, int, int, int, int);
/*
   This is an example of using scatterv to send a matrix from one
   process to all others, with the matrix stored in Fortran order.
   Note the use of an explicit UB to enable the sources to overlap.

   This tests scatterv to make sure that it uses the datatype size
   and extent correctly.  It requires number of processors that
   can be split with MPI_Dims_create.

 */

void SetData(double *sendbuf, double *recvbuf, int nx, int ny,
             int myrow, int mycol, int nrow, int ncol)
{
    int coldim, i, j, m, k;
    double *p;

    if (myrow == 0 && mycol == 0) {
        coldim = nx * nrow;
        for (j = 0; j < ncol; j++) {
            for (i = 0; i < nrow; i++) {
                p = sendbuf + i * nx + j * (ny * coldim);
                for (m = 0; m < ny; m++) {
                    for (k = 0; k < nx; k++) {
                        p[k] = 1000 * j + 100 * i + m * nx + k;
                    }
                    p += coldim;
                }
            }
        }
    }
    for (i = 0; i < nx * ny; i++)
        recvbuf[i] = -1.0;
}

int CheckData(double *recvbuf, int nx, int ny, int myrow, int mycol, int nrow, int expect_no_value)
{
    int coldim, m, k;
    double *p, val;
    int errs = 0;

    coldim = nx;
    p = recvbuf;
    for (m = 0; m < ny; m++) {
        for (k = 0; k < nx; k++) {
            /* If expect_no_value is true then we assume that the pre-scatterv
             * value should remain in the recvbuf for our portion of the array.
             * This is the case for the root process when using MPI_IN_PLACE. */
            if (expect_no_value)
                val = -1.0;
            else
                val = 1000 * mycol + 100 * myrow + m * nx + k;

            if (p[k] != val) {
                errs++;
                if (errs < 10) {
                    printf("Error in (%d,%d) [%d,%d] location, got %f expected %f\n",
                           m, k, myrow, mycol, p[k], val);
                } else if (errs == 10) {
                    printf("Too many errors; suppressing printing\n");
                }
            }
        }
        p += coldim;
    }
    return errs;
}

int main(int argc, char **argv)
{
    int rank, size, myrow, mycol, nx, ny, stride, cnt, i, j, errs, errs_in_place;
    double *sendbuf, *recvbuf;
    MPI_Datatype vec, block;
    int *scdispls;
    MPI_Comm comm2d;
    int dims[2], periods[2], coords[2], lcoords[2];
    int *sendcounts;


    MTest_Init(&argc, &argv);
    MPI_Comm_rank(MPI_COMM_WORLD, &rank);
    MPI_Comm_size(MPI_COMM_WORLD, &size);

    /* Get a 2-d decomposition of the processes */
    dims[0] = 0;
    dims[1] = 0;
    MPI_Dims_create(size, 2, dims);
    periods[0] = 0;
    periods[1] = 0;
    MPI_Cart_create(MPI_COMM_WORLD, 2, dims, periods, 0, &comm2d);
    MPI_Cart_get(comm2d, 2, dims, periods, coords);
    myrow = coords[0];
    mycol = coords[1];
/*
    if (rank == 0)
        printf("Decomposition is [%d x %d]\n", dims[0], dims[1]);
*/

    /* Get the size of the matrix */
    nx = 10;
    ny = 8;
    stride = nx * dims[0];

    recvbuf = (double *) malloc(nx * ny * sizeof(double));
    if (!recvbuf) {
        MPI_Abort(MPI_COMM_WORLD, 1);
    }
    sendbuf = 0;
    if (myrow == 0 && mycol == 0) {
        sendbuf = (double *) malloc(nx * ny * size * sizeof(double));
        if (!sendbuf) {
            MPI_Abort(MPI_COMM_WORLD, 1);
        }
    }
    sendcounts = (int *) malloc(size * sizeof(int));
    scdispls = (int *) malloc(size * sizeof(int));

    MPI_Type_vector(ny, nx, stride, MPI_DOUBLE, &vec);
    MPI_Type_create_resized(vec, 0, nx * sizeof(double), &block);
    MPI_Type_free(&vec);
    MPI_Type_commit(&block);

    /* Set up the transfer */
    cnt = 0;
    for (i = 0; i < dims[1]; i++) {
        for (j = 0; j < dims[0]; j++) {
            sendcounts[cnt] = 1;
            /* Using Cart_coords makes sure that ranks (used by
             * sendrecv) matches the cartesian coordinates (used to
             * set data in the matrix) */
            MPI_Cart_coords(comm2d, cnt, 2, lcoords);
            scdispls[cnt++] = lcoords[0] + lcoords[1] * (dims[0] * ny);
        }
    }

    SetData(sendbuf, recvbuf, nx, ny, myrow, mycol, dims[0], dims[1]);
    MPI_Scatterv(sendbuf, sendcounts, scdispls, block, recvbuf, nx * ny, MPI_DOUBLE, 0, comm2d);
    if ((errs = CheckData(recvbuf, nx, ny, myrow, mycol, dims[0], 0))) {
        fprintf(stdout, "Failed to transfer data\n");
    }

    /* once more, but this time passing MPI_IN_PLACE for the root */
    SetData(sendbuf, recvbuf, nx, ny, myrow, mycol, dims[0], dims[1]);
    MPI_Scatterv(sendbuf, sendcounts, scdispls, block,
                 (rank == 0 ? MPI_IN_PLACE : recvbuf), nx * ny, MPI_DOUBLE, 0, comm2d);
    errs_in_place = CheckData(recvbuf, nx, ny, myrow, mycol, dims[0], (rank == 0));
    if (errs_in_place) {
        fprintf(stdout, "Failed to transfer data (MPI_IN_PLACE)\n");
    }

    errs += errs_in_place;

    if (sendbuf)
        free(sendbuf);
    free(recvbuf);
    free(sendcounts);
    free(scdispls);
    MPI_Type_free(&block);
    MPI_Comm_free(&comm2d);
    MTest_Finalize(errs);
    return MTestReturnValue(errs);
}