SIONlib  1.7.1
Scalable I/O library for parallel access to task-local files
sion_ompi_internal_gen.c
Go to the documentation of this file.
1 /****************************************************************************
2 ** SIONLIB http://www.fz-juelich.de/jsc/sionlib **
3 *****************************************************************************
4 ** Copyright (c) 2008-2016 **
5 ** Forschungszentrum Juelich, Juelich Supercomputing Centre **
6 ** **
7 ** See the file COPYRIGHT in the package base directory for details **
8 ****************************************************************************/
9 
14 #include <stdlib.h>
15 #include <stdio.h>
16 #include <stdarg.h>
17 #include <string.h>
18 #include <time.h>
19 
20 #include "mpi.h"
21 #include "omp.h"
22 
23 #if defined(_BGL)
24 #include <rts.h>
25 #ifndef __USE_FILE_OFFSET64
26 #define __USE_FILE_OFFSET64
27 #endif
28 #endif
29 
30 #include <sys/time.h>
31 
32 #include <sys/types.h>
33 #include <fcntl.h>
34 
35 #include <unistd.h>
36 
37 #include "sion.h"
38 #include "sion_debug.h"
39 #include "sion_error_handler.h"
40 #include "sion_internal.h"
41 #include "sion_fd.h"
42 #include "sion_filedesc.h"
43 #include "sion_printts.h"
44 
45 #include "sion_ompi_internal_gen.h"
46 
47 #ifdef SION_OMPI
48 
49 #include "sion_ompi.h"
50 
51 int _sion_gen_info_from_gcomm_ompi(int numFiles, MPI_Comm gComm, int *filenumber, int *lrank, int *lsize)
52 {
53  int gtasks, gRank;
54  int rc = SION_SUCCESS;
55  int task_per_file;
56 
57 
58  MPI_Comm_size(gComm, &gtasks);
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));
62 
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",
66  numFiles, gtasks));
67  }
68 
69  task_per_file = gtasks / numFiles;
70 
71  /* remaining tasks are added to last communicator */
72  if (gRank >= (numFiles * task_per_file)) {
73  *filenumber = numFiles - 1;
74  }
75  else {
76  *filenumber = gRank / task_per_file;
77  }
78 
79  if(*filenumber == numFiles - 1) {
80  *lrank=gRank-(numFiles - 1)*task_per_file;
81  *lsize=gtasks-(numFiles - 1)*task_per_file;
82  } else {
83  *lrank=gRank%task_per_file;
84  *lsize=task_per_file;
85  }
86 
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));
89 
90  return (rc);
91 }
92 
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;
96  int ncomms, color;
97  int rc = SION_SUCCESS;
98  MPI_Status status;
99 
100  MPI_Comm_size(gComm, &gSize);
101  MPI_Comm_rank(gComm, &gRank);
102 
103  if (lComm != MPI_COMM_NULL) {
104  MPI_Comm_size(lComm, &lSize);
105  MPI_Comm_rank(lComm, &lRank);
106  }
107  else {
108  lSize = gSize;
109  lRank = gRank;
110  }
111 
112  DPRINTFP((32, "_sion_get_info_from_splitted_comm_mpi", gRank, "lRank: %d\n", lRank));
113 
114  if (gRank == 0) {
115 
116  allRanks = (int *) malloc(gSize * sizeof(int));
117  allSizes = (int *) malloc(gSize * sizeof(int));
118  if ((allRanks == NULL) || (allSizes == NULL)) {
119  free(allRanks);
120  free(allSizes);
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)));
122  }
123  }
124 
125  /* Gather all local ranks to task 0 for counting.. */
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);
128 
129  ncomms = 0;
130  ntasks = 0;
131  if (gRank == 0) {
132  for (i = 0; i < gSize; i++) {
133  if (allRanks[i] == 0) {
134  /* Use the current number of the communicator for the file suffix!! */
135  if (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));
138  }
139  else {
140  /* it's my self */
141  color = ncomms;
142  }
143  ncomms++;
144  ntasks += allSizes[i];
145  }
146  }
147  }
148  /* Recv the num of communicator from global 0 =>used for the suffix */
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));
152 
153  }
154 
155  MPI_Bcast(&ncomms, 1, MPI_INT, 0, gComm);
156  MPI_Bcast(&ntasks, 1, MPI_INT, 0, gComm);
157 
158 
159  DPRINTFP((32, "_sion_get_info_from_splitted_comm_mpi", gRank, "#Comms=%05d #Tasks=%05d #TotalTasks=%05d\n", ncomms, ntasks, gSize));
160 
161  if (lComm != MPI_COMM_NULL) {
162  MPI_Bcast(&color, 1, MPI_INT, 0, lComm);
163  }
164 
165  if (gRank == 0) {
166  if(allRanks) free(allRanks);
167  if(allSizes) free(allSizes);
168  }
169  /* return parameters */
170  *numComm = ncomms;
171  *CommNumber = color;
172  *lrank = lRank;
173  *lsize = lSize;
174 
175  return (rc);
176 }
177 
178 int _sion_errorprint_ompi(int rc, int level, const char *format, ...)
179 {
180  int rank=-1;
181  va_list ap;
182 
183  int thread = omp_get_thread_num();
184  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
185 
186  va_start(ap, format);
187  rc=__sion_errorprint_vargs(rc, level, rank, thread, format, ap);
188  va_end(ap);
189 
190  return (rc);
191 }
192 
193 
194 /* Help Functions */
195 
196 int _sion_map_rank_ompi_to_mpi(int ompi_rank, int num_threads) {
197  int mpi_rank = 0;
198 
199  mpi_rank = (int) (ompi_rank / num_threads);
200 
201  return mpi_rank;
202 }
203 
204 int _sion_map_rank_ompi_to_thread_num(int ompi_rank, int num_threads) {
205  int thread_num = 0;
206 
207  thread_num = (int) (ompi_rank % num_threads);
208 
209  return thread_num;
210 }
211 
212 
213 int _sion_map_rank_mpi_to_ompi(int mpi_rank, int num_threads, int thread_num) {
214  int ompi_rank = 0;
215 
216  ompi_rank = mpi_rank * num_threads + thread_num;
217 
218  return ompi_rank;
219 }
220 
221 int _sion_get_size_ompi(int mpi_size, int num_threads) {
222  int ompi_size = 0;
223 
224  ompi_size = mpi_size * num_threads;
225 
226  return ompi_size;
227 }
228 
229 
230 /* end of ifdef OMPI */
231 #endif
Sion Time Stamp Header.