24 #ifndef __USE_FILE_OFFSET64 25 #define __USE_FILE_OFFSET64 31 #include <sys/types.h> 38 #include "sion_error_handler.h" 50 int _sion_errorprint_mpi(
int rc,
int level,
const char *format, ...)
52 int rank=-1, thread=-1;
56 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
60 rc=__sion_errorprint_vargs(rc, level, rank, thread, format, ap);
82 int rc = SION_SUCCESS;
86 MPI_Comm_size(gComm, >asks);
87 MPI_Comm_rank(gComm, &gRank);
89 if (gtasks < numFiles) {
90 return(_sion_errorprint_mpi(SION_NOT_SUCCESS,_SION_ERROR_RETURN,
91 "_sion_gen_info_from_gcomm_mpi: Number of files(%d) cannot be bigger the the number of tasks(%d), aborting...\n",
95 task_per_file = gtasks / numFiles;
98 if (gRank >= (numFiles * task_per_file)) {
99 *filenumber = numFiles - 1;
102 *filenumber = gRank / task_per_file;
105 if(*filenumber == numFiles - 1) {
106 *lrank=gRank-(numFiles - 1)*task_per_file;
107 *lsize=gtasks-(numFiles - 1)*task_per_file;
109 *lrank=gRank%task_per_file;
110 *lsize=task_per_file;
113 DPRINTFP((1,
"_sion_get_info_from_splitted_comm_mpi", gRank,
"Global communicator divided in %d local communicators (%f: %d of %d)\n",
114 numFiles,*filenumber,*lrank,*lsize));
119 int _sion_get_info_from_splitted_comm_mpi(MPI_Comm gComm, MPI_Comm lComm,
int *numComm,
int *CommNumber,
int *lrank,
int *lsize) {
120 int gSize, gRank, lSize, lRank;
121 int *allRanks=NULL, *allSizes=NULL, i, ntasks;
122 int ncomms, color = 0;
123 int rc = SION_SUCCESS;
126 MPI_Comm_size(gComm, &gSize);
127 MPI_Comm_rank(gComm, &gRank);
129 if (lComm != MPI_COMM_NULL) {
130 MPI_Comm_size(lComm, &lSize);
131 MPI_Comm_rank(lComm, &lRank);
138 DPRINTFP((32,
"_sion_get_info_from_splitted_comm_mpi", gRank,
"lRank: %d\n", lRank));
142 allRanks = (
int *) malloc(gSize *
sizeof(
int));
143 allSizes = (
int *) malloc(gSize *
sizeof(
int));
144 if ((allRanks == NULL) || (allSizes == NULL)) {
147 return(_sion_errorprint_mpi(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)));
152 MPI_Gather(&lRank, 1, MPI_INT, allRanks, 1, MPI_INT, 0, gComm);
153 MPI_Gather(&lSize, 1, MPI_INT, allSizes, 1, MPI_INT, 0, gComm);
158 for (i = 0; i < gSize; i++) {
159 if (allRanks[i] == 0) {
162 MPI_Send(&ncomms, 1, MPI_INT, i, i, gComm);
163 DPRINTFP((32,
"sion_paropen_comms_mpi", gRank,
"Sent ncomms=%05d to %d with TAG %d\n", ncomms, i, i));
170 ntasks += allSizes[i];
175 if ((lRank == 0) && (gRank != 0)) {
176 MPI_Recv(&color, 1, MPI_INT, 0, gRank, gComm, &status);
177 DPRINTFP((32,
"sion_paropen_comms_mpi", gRank,
"Received ncomms=%05d from %d with TAG %d\n", color, status.MPI_SOURCE, status.MPI_TAG));
181 MPI_Bcast(&ncomms, 1, MPI_INT, 0, gComm);
182 MPI_Bcast(&ntasks, 1, MPI_INT, 0, gComm);
185 DPRINTFP((32,
"_sion_get_info_from_splitted_comm_mpi", gRank,
"#Comms=%05d #Tasks=%05d #TotalTasks=%05d\n", ncomms, ntasks, gSize));
187 if (lComm != MPI_COMM_NULL) {
188 MPI_Bcast(&color, 1, MPI_INT, 0, lComm);
192 if(allRanks) free(allRanks);
193 if(allSizes) free(allSizes);
204 int _sion_get_numfiles_from_file_mpi(
char *fname) {
209 sion_int64 chunksize;
210 sion_int32 fsblksize;
213 sid = _sion_open_rank(fname,
"br", &chunksize, &fsblksize, &rank, &fileptr);
215 numfiles = sion_get_number_of_files(sid);
217 _sion_close_sid(sid);
219 DPRINTFP((1,
"_sion_get_numfiles_from_file_mpi", _SION_DEFAULT_RANK,
"sion file %s has %d files\n", fname, numfiles));
224 int _sion_get_filenumber_from_files_mpi(
char *fname, MPI_Comm gComm,
int *numfiles,
int *filenumber,
int *lRank) {
226 int sid, gSize, gRank, ntasks, nfiles;
227 int rc = SION_SUCCESS;
229 sion_int64 *chunksizes;
230 sion_int32 fsblksize;
232 int mapping_size = 0;
233 sion_int32 *mapping = NULL;
236 MPI_Comm_size(gComm, &gSize);
237 MPI_Comm_rank(gComm, &gRank);
239 DPRINTFP((1,
"_sion_get_filenumber_from_files_mpi", gRank,
"before open\n"));
242 chunksizes=NULL;globalranks=NULL;
243 sid=
_sion_open_read(fname,_SION_FMODE_READ|_SION_FMODE_ANSI,_SION_READ_MASTER_ONLY_OF_MULTI_FILES,
244 &ntasks,&nfiles,&chunksizes,&fsblksize,&globalranks,&fileptr);
247 DPRINTFP((1,
"_sion_get_filenumber_from_files_mpi", gRank,
"after open\n"));
249 DPRINTFP((1,
"_sion_get_filenumber_from_files_mpi", gRank,
"sion file %d files rc=%d\n", *numfiles, rc));
256 MPI_Bcast(numfiles, 1, MPI_INT, 0, gComm);
258 if((gRank == 0) && (*numfiles>1)) {
259 if(mapping_size!=gSize) {
260 return(_sion_errorprint_mpi(SION_NOT_SUCCESS,_SION_ERROR_RETURN,
"_sion_get_filenumber_from_files_mpi: Incorrect sum of ntasks of files %d <> %d\n", mapping_size, gSize));
265 return(_sion_errorprint_mpi(SION_NOT_SUCCESS,_SION_ERROR_RETURN,
"_sion_get_filenumber_from_files_mpi: could not get numfiles from sion file\n"));
269 DPRINTFP((1,
"_sion_get_filenumber_from_files_mpi", gRank,
"before scatter\n"));
270 MPI_Scatter(mapping, 2, SION_MPI_INT32, lpos, 2, SION_MPI_INT32, 0, gComm);
273 DPRINTFP((1,
"_sion_get_filenumber_from_files_mpi", gRank,
"after scatter\n"));
277 DPRINTFP((1,
"_sion_get_filenumber_from_files_mpi", gRank,
"only one file -> filenumber=%d lRank=%d\n",*filenumber,*lRank));
282 if (sid>=0) _sion_close_sid(sid);
285 DPRINTFP((1,
"_sion_get_filenumber_from_files_mpi", gRank,
"close rc=%d\n",rc));
int sion_get_mapping(int sid, int *mapping_size, sion_int32 **mapping, int *numfiles)
Returns pointers to the internal field mapping.
int _sion_gen_info_from_gcomm_mpi(int numFiles, MPI_Comm gComm, int *filenumber, int *lrank, int *lsize)
Splits a Communicator in numfiles different communicators.
int _sion_open_read(const char *fname, sion_int64 file_mode_flags, int read_all, int *ntasks, int *nfiles, sion_int64 **chunksizes, sion_int32 *fsblksize, int **globalranks, FILE **fileptr)
internal sion serial open function for reading on one or more files