Skip to content
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion mpi-proxy-split/record-replay.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -306,7 +306,7 @@ restoreTypeCreateStruct(MpiRecord& rec)
int count = rec.args(0);
int *blocklengths = rec.args(1);
MPI_Aint *displs = rec.args(2);
MPI_Datatype *types = (MPI_Datatype*)(intptr_t)rec.args(3);
MPI_Datatype *types = (MPI_Datatype*)rec.args(3);
MPI_Datatype newtype = MPI_DATATYPE_NULL;
retval = FNC_CALL(Type_create_struct, rec)(count, blocklengths,
displs, types, &newtype);
Expand Down
2 changes: 1 addition & 1 deletion mpi-proxy-split/record-replay.h
Original file line number Diff line number Diff line change
Expand Up @@ -437,7 +437,7 @@ namespace dmtcp_mpi
{
MPI_Datatype newtype = (MPI_Datatype)(int)rec->args(4);
int count = rec->args(0);
MPI_Datatype *oldtypes = (MPI_Datatype*)(intptr_t)rec->args(3);
MPI_Datatype *oldtypes = (MPI_Datatype*)rec->args(3);
datatype_create(newtype);
datatype_incRef(count, oldtypes);
break;
Expand Down
2 changes: 1 addition & 1 deletion mpi-proxy-split/test/Makefile
Original file line number Diff line number Diff line change
Expand Up @@ -40,7 +40,7 @@ FILES=mpi_hello_world \
Abort_test Allreduce_test Alltoall_test Alltoallv_test \
Allgather_test Group_size_rank Type_commit_contiguous \
Irecv_test Alloc_mem sendrecv_replace_test \
Type_hvector_test Type_vector_test \
Type_hvector_test Type_vector_test Type_create_struct_test \
${FILES_FORTRAN} \
send_recv send_recv_loop large_async_p2p \
File_read_write_test File_characteristics_test File_size_test \
Expand Down
83 changes: 83 additions & 0 deletions mpi-proxy-split/test/Type_create_struct_test.c
Original file line number Diff line number Diff line change
@@ -0,0 +1,83 @@
/*
Test for the MPI_Type_create_struct method

Run with any ranks
Defaults to 10000 iterations
Intended to be run with mana_test.py

Inspired by:
http://mpi.deino.net/mpi_functions/MPI_Type_create_struct.html
*/

#include <assert.h>
#include <mpi.h>
#include <stdio.h>
#include <stdlib.h>
#include <stddef.h>
#include <string.h>

struct Partstruct {
char c;
double d[6];
char b[7];
};

int main(int argc, char **argv) {
int max_iterations = 10000;
if (argc > 1) max_iterations = atoi(argv[1]);

int rank, comm_size;
MPI_Init(&argc, &argv);
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
MPI_Comm_size(MPI_COMM_WORLD, &comm_size);

if (rank == 0) {
printf("Running Partstruct test with %d ranks and %d iterations\n", comm_size, max_iterations);
}

MPI_Datatype partstruct_type;

int block_lengths[3] = {1, 6, 7};
MPI_Datatype types[3] = {MPI_CHAR, MPI_DOUBLE, MPI_CHAR};
MPI_Aint offsets[3];

offsets[0] = offsetof(struct Partstruct, c);
offsets[1] = offsetof(struct Partstruct, d);
offsets[2] = offsetof(struct Partstruct, b);

MPI_Type_create_struct(3, block_lengths, offsets, types, &partstruct_type);
MPI_Type_commit(&partstruct_type);

for (int i = 0; i < max_iterations; i++) {
struct Partstruct send_p;
struct Partstruct recv_p;

send_p.c = 'A' + (i % 26);
for (int j = 0; j < 6; j++) send_p.d[j] = (double)i + j;
snprintf(send_p.b, 7, "iter%d", i % 100);

int dest = (rank + 1) % comm_size;
int source = (rank - 1 + comm_size) % comm_size;

MPI_Sendrecv(&send_p, 1, partstruct_type, dest, 0,
&recv_p, 1, partstruct_type, source, 0,
MPI_COMM_WORLD, MPI_STATUS_IGNORE);

if (comm_size == 1) {
assert(recv_p.c == send_p.c);
for (int j = 0; j < 6; j++) assert(recv_p.d[j] == send_p.d[j]);
assert(strcmp(recv_p.b, send_p.b) == 0);
}

if (rank == 0 && i % 1000 == 999) {
printf("Completed %d iterations\n", i + 1);
fflush(stdout);
}
}

if (rank == 0) printf("Test finished successfully.\n");

MPI_Type_free(&partstruct_type);
MPI_Finalize();
return EXIT_SUCCESS;
}