14 #define _XOPEN_SOURCE 700
26 #define MPI_Comm_rank PMPI_Comm_rank
27 #define MPI_Comm_size PMPI_Comm_size
28 #define MPI_Gather PMPI_Gather
29 #define MPI_Scatter PMPI_Scatter
30 #define MPI_Bcast PMPI_Bcast
31 #define MPI_Barrier PMPI_Barrier
32 #define MPI_Comm_split PMPI_Comm_split
33 #define MPI_Send PMPI_Send
34 #define MPI_Recv PMPI_Recv
39 #include <sys/types.h>
46 #include "sion_error_handler.h"
56 int _sion_register_callbacks_mpi(
void) {
58 aid=sion_generic_create_api(
"SIONlib_MPI_API");
61 sion_generic_register_create_local_commgroup_cb(aid,&_sion_mpi_create_lcg_cb);
62 sion_generic_register_free_local_commgroup_cb(aid,&_sion_mpi_free_lcg_cb);
64 sion_generic_register_barrier_cb(aid,&_sion_mpi_barrier_cb);
65 sion_generic_register_bcastr_cb(aid,&_sion_mpi_bcastr_cb);
66 sion_generic_register_gatherr_cb(aid,&_sion_mpi_gatherr_cb);
67 sion_generic_register_scatterr_cb(aid,&_sion_mpi_scatterr_cb);
68 sion_generic_register_gathervr_cb(aid,&_sion_mpi_gathervr_cb);
69 sion_generic_register_scattervr_cb(aid,&_sion_mpi_scattervr_cb);
70 sion_generic_register_gather_and_execute_cb(aid,&_sion_mpi_gather_process_cb);
71 sion_generic_register_execute_and_scatter_cb(aid,&_sion_mpi_process_scatter_cb);
72 sion_generic_register_get_capability_cb(aid,&_sion_mpi_get_capability_cb);
77 int _sion_mpi_create_lcg_cb(
void **local_commgroup,
void *global_commgroup,
80 int filenumber,
int numfiles
86 int create_lcomm=1, set_in_global=1, mrank=0, msize=1, color;
88 DPRINTFP((256,
"_mpi_create_lcg_cb", grank,
" split now comm: grank=%d gsize=%d filenumber=%d, numfiles=%d, lrank=%d lsize=%d \n",
89 grank, gsize, filenumber, numfiles, lrank, lsize));
91 if(global_commgroup==NULL) {
92 fprintf(stderr,
"_mpi_create_lcg_cb: error no global commgroup given, aborting ...\n");
95 if(*local_commgroup!=NULL) {
96 fprintf(stderr,
"_mpi_create_lcg_cb: error local commgroup already initialized, aborting ...\n");
101 if(sapi_global->lcommgroup!=NULL) {
103 if(sapi_global->lcommgroup->commset==0) {
104 *local_commgroup=sapi_global->lcommgroup;
105 create_lcomm=0;set_in_global=0;
106 sapi_global->lcommgroup->commset=1;
108 create_lcomm=1;set_in_global=0;
116 if (commgroup == NULL) {
117 fprintf(stderr,
"_mpi_create_lcg_cb: cannot allocate memory for local commgroup of size %lu, aborting ...\n",
122 if(filenumber==-1) color=MPI_UNDEFINED;
124 rc = MPI_Comm_split(sapi_global->comm, color, lrank, &commgroup->comm);
125 DPRINTFP((256,
"_mpi_create_lcg_cb", grank,
" rc=%d from MPI_Comm_split(comm,%d,%d,&newcomm)\n",rc,color,lrank));
126 commgroup->local=1; commgroup->commset=1; commgroup->lcommgroup=NULL;
127 commgroup->commcreated=1;
128 commgroup->rank=lrank;
129 commgroup->size=lsize;
134 sapi_global->lcommgroup=commgroup;
137 *local_commgroup=commgroup;
139 if((filenumber!=-1) && commgroup) {
140 MPI_Comm_rank(commgroup->comm, &mrank);
141 MPI_Comm_size(commgroup->comm, &msize);
144 DPRINTFP((256,
"_mpi_create_lcg_cb", grank,
" leave rc=%d rank %d of %d\n",rc,mrank,msize));
149 int _sion_mpi_free_lcg_cb(
void *local_commgroup) {
153 if ( (commgroup->commset) && (commgroup->commcreated) ) {
154 DPRINTFP((256,
"_mpi_free_lcg_cb", commgroup->rank,
" free now comm\n"));
155 rc=MPI_Comm_free(&commgroup->comm);
156 DPRINTFP((256,
"_mpi_free_lcg_cb", commgroup->rank,
" free now comm rc=%d\n",rc));
162 int _sion_mpi_barrier_cb(
void *commdata)
166 MPI_Comm commp = sapi->comm;
167 rc = MPI_Barrier(commp);
171 int _sion_mpi_bcastr_cb(
void *data,
void *commdata,
int dtype,
int nelem,
int root)
175 MPI_Comm commp = sapi->comm;
178 rc = MPI_Bcast((sion_int32 *) data, nelem, SION_MPI_INT32, root, commp);
181 rc = MPI_Bcast((sion_int64 *) data, nelem, SION_MPI_INT64, root, commp);
184 rc = MPI_Bcast((
char *) data, nelem, MPI_CHAR, root, commp);
187 rc = MPI_Bcast((sion_int64 *) data, nelem, SION_MPI_INT64, root, commp);
193 int _sion_mpi_gatherr_cb(
void *indata,
void *outdata,
void *commdata,
int dtype,
int nelem,
int root)
198 MPI_Comm commp = sapi->comm;
200 MPI_Comm_rank(commp, &rank);
201 MPI_Comm_size(commp, &size);
203 DPRINTFP((256,
"_mpi_gatherr_cb", rank,
" gatherr on %d of %d nelem=%d root=%d\n", rank, size, nelem, root));
208 rc = MPI_Gather((sion_int32 *) indata, nelem, SION_MPI_INT32, (sion_int32 *) outdata, nelem, SION_MPI_INT32, root, commp);
211 rc = MPI_Gather((sion_int64 *) indata, nelem, SION_MPI_INT64, (sion_int64 *) outdata, nelem, SION_MPI_INT64, root, commp);
214 rc = MPI_Gather((
char *) indata, nelem, MPI_CHAR, (
char *) outdata, nelem, MPI_CHAR, root, commp);
217 rc = MPI_Gather((sion_int64 *) indata, nelem, SION_MPI_INT64, (sion_int64 *) outdata, nelem, SION_MPI_INT64, root, commp);
225 int _sion_mpi_scatterr_cb(
void *indata,
void *outdata,
void *commdata,
int dtype,
int nelem,
int root)
229 MPI_Comm commp = sapi->comm;
230 ONLY_DEBUG(
int rank=sapi->rank);
232 DPRINTFP((256,
"_mpi_scatterr_cb", rank,
" starting nelem=%d root=%d\n", nelem, root));
236 rc = MPI_Scatter((sion_int32 *) indata, nelem, SION_MPI_INT32, (sion_int32 *) outdata, nelem, SION_MPI_INT32, root, commp);
239 rc = MPI_Scatter((sion_int64 *) indata, nelem, SION_MPI_INT64, (sion_int64 *) outdata, nelem, SION_MPI_INT64, root, commp);
242 rc = MPI_Scatter((
char *) indata, nelem, MPI_CHAR, (
char *) outdata, nelem, MPI_CHAR, root, commp);
245 rc = MPI_Scatter((sion_int64 *) indata, nelem, SION_MPI_INT64, (sion_int64 *) outdata, nelem, SION_MPI_INT64, root, commp);
252 int _sion_mpi_gathervr_cb(
void *indata,
void *outdata,
void *commdata,
int dtype,
int *counts,
int nelem,
int root)
255 int size, rank, t, offset;
258 MPI_Comm commp = sapi->comm;
261 MPI_Comm_rank(commp, &rank);
262 MPI_Comm_size(commp, &size);
264 DPRINTFP((256,
"_mpi_gathervr_cb", rank,
" input nelem=%d root=%d indata=%x, outdata=%x\n", nelem, root, indata, outdata));
268 displs = (
int *) malloc(size *
sizeof(
int));
269 if (displs == NULL) {
270 fprintf(stderr,
"_mpi_gathervr_cb: cannot allocate temporary memory of size %zu (displs), aborting ...\n",(
size_t) size *
sizeof(
int));
274 for(t=0;t<size;t++) {
284 rc = MPI_Gatherv((sion_int32 *) indata, nelem, SION_MPI_INT32, (sion_int32 *) outdata, counts, displs, SION_MPI_INT32, root, commp);
287 rc = MPI_Gatherv((sion_int64 *) indata, nelem, SION_MPI_INT64, (sion_int64 *) outdata, counts, displs, SION_MPI_INT64, root, commp);
290 rc = MPI_Gatherv((
char *) indata, nelem, MPI_CHAR, (sion_int32 *) outdata, counts, displs, MPI_CHAR, root, commp);
293 rc = MPI_Gatherv((sion_int64 *) indata, nelem, SION_MPI_INT64, (sion_int64 *) outdata, counts, displs, SION_MPI_INT64, root, commp);
298 if(displs) free(displs);
303 int _sion_mpi_scattervr_cb(
void *indata,
void *outdata,
void *commdata,
int dtype,
int *counts,
int nelem,
int root)
306 int size, rank, t, offset;
309 MPI_Comm commp = sapi->comm;
312 MPI_Comm_rank(commp, &rank);
313 MPI_Comm_size(commp, &size);
315 DPRINTFP((256,
"_mpi_scattervr_cb", rank,
" input nelem=%d root=%d\n", nelem, root));
319 displs = (
int *) malloc(size *
sizeof(
int));
320 if (displs == NULL) {
321 fprintf(stderr,
"_mpi_scattervr_cb: cannot allocate temporary memory of size %zu (displs), aborting ...\n",(
size_t) size *
sizeof(
int));
325 for(t=0;t<size;t++) {
328 DPRINTFP((256,
"_mpi_scattervr_cb", rank,
" after MPI_Gather %2d -> dpls=%2d sendcounts=%d\n", t,displs[t],counts[t]));
335 rc = MPI_Scatterv((sion_int32 *) outdata, counts, displs, SION_MPI_INT32, (sion_int32 *) indata, nelem, SION_MPI_INT32, root, commp);
338 rc = MPI_Scatterv((sion_int64 *) outdata, counts, displs, SION_MPI_INT64, (sion_int64 *) indata, nelem, SION_MPI_INT64, root, commp);
341 rc = MPI_Scatterv((
char *) outdata, counts, displs, MPI_CHAR, (sion_int32 *) indata, nelem, MPI_CHAR, root, commp);
344 rc = MPI_Scatterv((sion_int64 *) outdata, counts, displs, SION_MPI_INT64, (sion_int64 *) indata, nelem, SION_MPI_INT64, root, commp);
349 if(displs) free(displs);