25 #ifndef __USE_FILE_OFFSET64 26 #define __USE_FILE_OFFSET64 32 #include <sys/types.h> 39 #include "sion_error_handler.h" 51 int _sion_gen_info_from_gcomm_ompi(
int numFiles, MPI_Comm gComm,
int *filenumber,
int *lrank,
int *lsize)
54 int rc = SION_SUCCESS;
58 MPI_Comm_size(gComm, >asks);
59 MPI_Comm_rank(gComm, &gRank);
60 DPRINTFP((1,
"_sion_get_info_from_splitted_comm_mpi", gRank,
"enter: gcomm: %d of %d, numfiles=%d\n",
61 gRank,gtasks,numFiles));
63 if (gtasks < numFiles) {
64 return(_sion_errorprint_ompi(SION_NOT_SUCCESS,_SION_ERROR_RETURN,
65 "_sion_gen_info_from_gcomm_mpi: Number of files(%d) cannot be bigger the the number of tasks(%d), aborting...\n",
69 task_per_file = gtasks / numFiles;
72 if (gRank >= (numFiles * task_per_file)) {
73 *filenumber = numFiles - 1;
76 *filenumber = gRank / task_per_file;
79 if(*filenumber == numFiles - 1) {
80 *lrank=gRank-(numFiles - 1)*task_per_file;
81 *lsize=gtasks-(numFiles - 1)*task_per_file;
83 *lrank=gRank%task_per_file;
87 DPRINTFP((1,
"_sion_get_info_from_splitted_comm_mpi", gRank,
"Global communicator divided in %d local communicators (%d: %d of %d)\n",
88 numFiles,*filenumber,*lrank,*lsize));
93 int _sion_get_info_from_splitted_comm_ompi(MPI_Comm gComm, MPI_Comm lComm,
int *numComm,
int *CommNumber,
int *lrank,
int *lsize) {
94 int gSize, gRank, lSize, lRank;
95 int *allRanks=NULL, *allSizes=NULL, i, ntasks;
97 int rc = SION_SUCCESS;
100 MPI_Comm_size(gComm, &gSize);
101 MPI_Comm_rank(gComm, &gRank);
103 if (lComm != MPI_COMM_NULL) {
104 MPI_Comm_size(lComm, &lSize);
105 MPI_Comm_rank(lComm, &lRank);
112 DPRINTFP((32,
"_sion_get_info_from_splitted_comm_mpi", gRank,
"lRank: %d\n", lRank));
116 allRanks = (
int *) malloc(gSize *
sizeof(
int));
117 allSizes = (
int *) malloc(gSize *
sizeof(
int));
118 if ((allRanks == NULL) || (allSizes == NULL)) {
121 return(_sion_errorprint_ompi(SION_NOT_SUCCESS,_SION_ERROR_RETURN,
"_sion_get_info_from_splitted_comm_mpi: Cannot allocate temp arrays of size %lu\n", (
unsigned long) gSize *
sizeof(
int)));
126 MPI_Gather(&lRank, 1, MPI_INT, allRanks, 1, MPI_INT, 0, gComm);
127 MPI_Gather(&lSize, 1, MPI_INT, allSizes, 1, MPI_INT, 0, gComm);
132 for (i = 0; i < gSize; i++) {
133 if (allRanks[i] == 0) {
136 MPI_Send(&ncomms, 1, MPI_INT, i, i, gComm);
137 DPRINTFP((32,
"sion_paropen_comms_mpi", gRank,
"Sent ncomms=%05d to %d with TAG %d\n", ncomms, i, i));
144 ntasks += allSizes[i];
149 if ((lRank == 0) && (gRank != 0)) {
150 MPI_Recv(&color, 1, MPI_INT, 0, gRank, gComm, &status);
151 DPRINTFP((32,
"sion_paropen_comms_mpi", gRank,
"Received ncomms=%05d from %d with TAG %d\n", color, status.MPI_SOURCE, status.MPI_TAG));
155 MPI_Bcast(&ncomms, 1, MPI_INT, 0, gComm);
156 MPI_Bcast(&ntasks, 1, MPI_INT, 0, gComm);
159 DPRINTFP((32,
"_sion_get_info_from_splitted_comm_mpi", gRank,
"#Comms=%05d #Tasks=%05d #TotalTasks=%05d\n", ncomms, ntasks, gSize));
161 if (lComm != MPI_COMM_NULL) {
162 MPI_Bcast(&color, 1, MPI_INT, 0, lComm);
166 if(allRanks) free(allRanks);
167 if(allSizes) free(allSizes);
178 int _sion_errorprint_ompi(
int rc,
int level,
const char *format, ...)
183 int thread = omp_get_thread_num();
184 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
186 va_start(ap, format);
187 rc=__sion_errorprint_vargs(rc, level, rank, thread, format, ap);
196 int _sion_map_rank_ompi_to_mpi(
int ompi_rank,
int num_threads) {
199 mpi_rank = (int) (ompi_rank / num_threads);
204 int _sion_map_rank_ompi_to_thread_num(
int ompi_rank,
int num_threads) {
207 thread_num = (int) (ompi_rank % num_threads);
213 int _sion_map_rank_mpi_to_ompi(
int mpi_rank,
int num_threads,
int thread_num) {
216 ompi_rank = mpi_rank * num_threads + thread_num;
221 int _sion_get_size_ompi(
int mpi_size,
int num_threads) {
224 ompi_size = mpi_size * num_threads;