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
|
/*
* Copyright (C) by Argonne National Laboratory
* See COPYRIGHT in top-level directory
*/
/*
* This program was sent in as an example that did not perform as expected.
* The program has a bug in that it is sending 3 characters but receiving
* three integers, which is not a valid MPI program (the type signatures
* must match). However, a good MPI implementation will handle this
* gracefully, which is why this test is included in the error directory
*/
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include "mpi.h"
#include "mpitest.h"
#define BUFSIZE 128
#define FAIL_ERROR 123
#define RESFAIL_ERROR 123
#define INTERNAL_ERROR 123
void if_error(const char *function, const char *data, int ret);
void if_error(const char *function, const char *data, int ret)
{
if (ret == 0)
return;
if (data)
printf("%s for %s returned %d (%#x)\n", function, data, ret, ret);
else
printf("%s returned %d (%#x)\n", function, ret, ret);
exit(INTERNAL_ERROR);
}
int main(int argc, char *argv[])
{
int ret, errs = 0;
char *src, *sendrec;
int bufsize = BUFSIZE;
int myrank, nprocs;
int i;
MPI_Status status;
int small_non_contig_struct_count = 3;
int small_non_contig_struct_blocklens[] = { 1, 1, 1 };
MPI_Aint small_non_contig_struct_disps[] = { 0, 2, 4 };
MPI_Datatype small_non_contig_struct_types[] = { MPI_CHAR, MPI_CHAR, MPI_CHAR };
MPI_Datatype small_non_contig_struct_type;
int contig_indexed_count = 3;
int contig_indexed_blocklens[] = { 1, 2, 1 };
int contig_indexed_indices[] = { 4, 8, 16 };
MPI_Datatype contig_indexed_inner_type = MPI_INT;
MPI_Datatype contig_indexed_type;
MTest_Init(&argc, &argv);
ret = MPI_Comm_rank(MPI_COMM_WORLD, &myrank);
ret = MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
if (nprocs < 2) {
printf("Need at least 2 procs\n");
exit(RESFAIL_ERROR);
}
ret = MPI_Type_create_struct(small_non_contig_struct_count,
small_non_contig_struct_blocklens,
small_non_contig_struct_disps,
small_non_contig_struct_types, &small_non_contig_struct_type);
if_error("MPI_Type_create_struct", "small_non_contig_struct_type", ret);
ret = MPI_Type_commit(&small_non_contig_struct_type);
if_error("MPI_Type_commit", "small_non_contig_struct_type", ret);
ret = MPI_Type_indexed(contig_indexed_count, contig_indexed_blocklens,
contig_indexed_indices, contig_indexed_inner_type, &contig_indexed_type);
if_error("MPI_Type_indexed", "contig_indexed_type", ret);
ret = MPI_Type_commit(&contig_indexed_type);
if_error("MPI_Type_commit", "contig_indexed_type", ret);
ret = MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &src);
if (ret != 0) {
printf("MPI_Alloc_mem src = #%x\n", ret);
exit(INTERNAL_ERROR);
}
ret = MPI_Alloc_mem(bufsize, MPI_INFO_NULL, &sendrec);
if (ret != 0) {
printf("MPI_Alloc_mem sendrec buf = #%x\n", ret);
exit(INTERNAL_ERROR);
}
for (i = 0; i < bufsize; i++) {
src[i] = (char) i + 1;
}
memset(sendrec, 0, bufsize);
MPI_Barrier(MPI_COMM_WORLD);
if (myrank == 1) {
MPI_Send(src, 1, small_non_contig_struct_type, 0, 0xabc, MPI_COMM_WORLD);
} else {
MPI_Comm_set_errhandler(MPI_COMM_WORLD, MPI_ERRORS_RETURN);
ret = MPI_Recv(sendrec, 1, contig_indexed_type, 1, 0xabc, MPI_COMM_WORLD, &status);
if (ret == MPI_SUCCESS) {
printf("MPI_Recv succeeded with non-matching datatype signature\n");
errs++;
}
}
MPI_Barrier(MPI_COMM_WORLD);
MPI_Type_free(&small_non_contig_struct_type);
MPI_Type_free(&contig_indexed_type);
MPI_Free_mem(src);
MPI_Free_mem(sendrec);
MTest_Finalize(errs);
return MTestReturnValue(errs);
}
|