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
|
/*
* Copyright (C) by Argonne National Laboratory
* See COPYRIGHT in top-level directory
*/
#include <stdio.h>
#include <string.h>
#include "mpi.h"
#include "mpitest.h"
#define X 64
#define Y 8
#define Z 512
double array[X][Y][Z];
int main(int argc, char *argv[])
{
int myrank;
MPI_Datatype subarray;
int array_size[] = { X, Y, Z };
int array_subsize[] = { X / 2, Y / 2, Z };
int array_start[] = { 0, 0, 0 };
int i, j, k;
int errs = 0;
MTest_Init(&argc, &argv);
MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
for (i = 0; i < X; ++i) {
for (j = 0; j < Y; ++j) {
for (k = 0; k < Z; ++k) {
if (myrank == 0)
array[i][j][k] = 2.0;
else
array[i][j][k] = -2.0;
}
}
}
MPI_Type_create_subarray(3, array_size, array_subsize, array_start, MPI_ORDER_C,
MPI_DOUBLE, &subarray);
MPI_Type_commit(&subarray);
if (myrank == 0)
MPI_Send(array, 1, subarray, 1, 0, MPI_COMM_WORLD);
else {
MPI_Recv(array, 1, subarray, 0, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
for (i = array_start[0]; i < array_subsize[0]; ++i) {
for (j = array_start[1]; j < array_subsize[1]; ++j) {
for (k = array_start[2]; k < array_subsize[2]; ++k) {
if (array[i][j][k] != 2.0)
++errs;
}
}
}
}
MPI_Type_free(&subarray);
MTest_Finalize(errs);
return MTestReturnValue(errs);
}
|