Nek5000
SEM for Incompressible NS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
nek_comm.c
Go to the documentation of this file.
1 /*
2  * MPI wrappers
3  *
4  * timer:
5  * 0 MPI_Allreduce
6  * 1 MPI_Allreduce(sync)
7  * 2 MPI_Waitall
8  * 3 MPI_Barrier
9  * 4 MPI_Irecv
10  * 5 MPI_Isend
11  * 6 MPI_Recv
12  * 7 MPI_Send
13  *
14  */
15 #include <stdio.h>
16 #ifdef MPI
17 #include <mpi.h>
18 #endif
19 #include <math.h>
20 
21 #define NTIMER 8
22 #define NCOUNTER NTIMER
23 
24 #ifdef UPCASE
25 # define FORTRAN_NAME(low,up) up
26 #else
27 #ifdef UNDERSCORE
28 # define FORTRAN_NAME(low,up) low##_
29 #else
30 # define FORTRAN_NAME(low,up) low
31 #endif
32 #endif
33 
34 #define nek_comm_settings FORTRAN_NAME(nek_comm_settings, NEK_COMM_SETTINGS)
35 #define nek_comm_getstat FORTRAN_NAME(nek_comm_getstat, NEK_COMM_GETSTAT)
36 #define nek_comm_startstat FORTRAN_NAME(nek_comm_startstat, NEK_COMM_STARTSTAT)
37 
38 
39 int SYNC = 0;
40 int TIMING = 1;
41 double TIMER[NTIMER] = {(double)0.0};
42 int COUNTER[NCOUNTER] = {0};
43 
44 
45 #if defined(MPITIMER)
46 
47 void nek_comm_settings(int *sync,int *timing)
48 {
49  SYNC = *sync;
50  TIMING = *timing;
51 }
52 
53 void nek_comm_getstat(double *timer, int *counter)
54 {
55  int i;
56  double t0,t1;
57 
58  /* estimate timer overhead */
59  t0 = MPI_Wtime();
60  t1 = MPI_Wtime();
61  while (t1 == t0) t1 = MPI_Wtime();
62  t1 -= t0;
63 
64 
65  /* substract timer overhead */
66 
67  /*
68  TIMER[0] -= 3*COUNTER[0] * t1;
69  TIMER[1] -= 2*COUNTER[1] * t1;
70  TIMER[2] -= 2*COUNTER[2] * t1;
71  TIMER[3] -= 2*COUNTER[3] * t1;
72  TIMER[4] -= 2*COUNTER[4] * t1;
73  TIMER[5] -= 2*COUNTER[5] * t1;
74  TIMER[6] -= 2*COUNTER[6] * t1;
75  TIMER[7] -= 2*COUNTER[7] * t1;
76  */
77 
78  TIMER[1] += TIMER[0];
79  COUNTER[1] = COUNTER[0];
80  for (i = 0; i < NTIMER; i++) timer[i] = fmax(TIMER[i],(double)0.0);
81  for (i = 0; i < NCOUNTER; i++) counter[i] = COUNTER[i];
82 }
83 
84 void nek_comm_startstat(void)
85 {
86  int i;
87  for (i = 0; i < NTIMER; i++) TIMER[i] = 0.0;
88  for (i = 0; i < NCOUNTER; i++) COUNTER[i] = 0;
89 }
90 
91 /* FORTRAN wrappers */
92 
93 
94 #pragma weak MPI_ALLREDUCE = mpi_allreduce_f
95 #pragma weak mpi_allreduce = mpi_allreduce_f
96 #pragma weak mpi_allreduce_ = mpi_allreduce_f
97 #pragma weak mpi_allreduce__ = mpi_allreduce_f
98 void mpi_allreduce_f(char *sendbuf, char *recvbuf, MPI_Fint *count,
100  MPI_Fint *ierr)
101 {
102  MPI_Comm c_comm = MPI_Comm_f2c(*comm);
103  MPI_Datatype c_type = MPI_Type_f2c(*datatype);
104  MPI_Op c_op = MPI_Op_f2c(*op);
105  double t0,t1;
106 
107  COUNTER[0]++;
108 
109  if (TIMING) t0 = MPI_Wtime();
110  if (SYNC) *ierr = MPI_Barrier(c_comm);
111  if (TIMING) t1 = MPI_Wtime();
112  *ierr = PMPI_Allreduce(sendbuf, recvbuf, *count, c_type, c_op, c_comm);
113  if (TIMING) {TIMER[0] += MPI_Wtime()-t1; TIMER[1] += t1-t0;}
114 }
115 
116 #pragma weak MPI_BARRIER = mpi_barrier_f
117 #pragma weak mpi_barrier = mpi_barrier_f
118 #pragma weak mpi_barrier_ = mpi_barrier_f
119 #pragma weak mpi_barrier__ = mpi_barrier_f
120 void mpi_barrier_f(MPI_Fint *comm, MPI_Fint *ierr)
121 {
122  MPI_Comm c_comm = MPI_Comm_f2c(*comm);
123  double t0;
124 
125  COUNTER[3]++;
126 
127  if (TIMING) t0 = MPI_Wtime();
128  *ierr = PMPI_Barrier(c_comm);
129  if (TIMING) TIMER[3] += MPI_Wtime()-t0;
130 }
131 
132 #pragma weak MPI_RECV = mpi_recv_f
133 #pragma weak mpi_recv = mpi_recv_f
134 #pragma weak mpi_recv_ = mpi_recv_f
135 #pragma weak mpi_recv__ = mpi_recv_f
136 void mpi_recv_f(char *buf, MPI_Fint *count, MPI_Fint *datatype,
137  MPI_Fint *source, MPI_Fint *tag, MPI_Fint *comm,
138  MPI_Fint *status, MPI_Fint *ierr)
139 {
140  MPI_Comm c_comm = MPI_Comm_f2c(*comm);
141  MPI_Datatype c_type = MPI_Type_f2c(*datatype);
142  MPI_Status *c_status = (MPI_Status *) status; /* FORTRAN_INTEGER == SIZEOF_INT */
143  double t0;
144 
145  COUNTER[6]++;
146 
147  if (TIMING) t0 = MPI_Wtime();
148  *ierr = PMPI_Recv(buf, *count, c_type, *source, *tag, c_comm, c_status);
149  if (TIMING) TIMER[6] += MPI_Wtime()-t0;
150 }
151 
152 #pragma weak MPI_SEND = mpi_send_f
153 #pragma weak mpi_send = mpi_send_f
154 #pragma weak mpi_send_ = mpi_send_f
155 #pragma weak mpi_send__ = mpi_send_f
156 void mpi_send_f(char *buf, MPI_Fint *count, MPI_Fint *datatype,
157  MPI_Fint *dest, MPI_Fint *tag, MPI_Fint *comm, MPI_Fint *ierr)
158 {
159  MPI_Comm c_comm = MPI_Comm_f2c(*comm);
160  MPI_Datatype c_type = MPI_Type_f2c(*datatype);
161  double t0;
162 
163  COUNTER[7]++;
164 
165  if (TIMING) t0 = MPI_Wtime();
166  *ierr = PMPI_Send(buf, *count, c_type, *dest, *tag, c_comm);
167  if (TIMING) TIMER[7] += MPI_Wtime()-t0;
168 }
169 
170 
171 /* C wrappers */
172 
173 
174 #pragma weak MPI_Allreduce = mpi_allreduce_c
175 #pragma weak MPI_Allreduce_ = mpi_allreduce_c
176 int mpi_allreduce_c(void *sendbuf, void *recvbuf, int count,
177  MPI_Datatype datatype, MPI_Op op, MPI_Comm comm)
178 {
179  int ierr;
180  double t0,t1;
181 
182  COUNTER[0]++;
183 
184  if (TIMING) t0 = MPI_Wtime();
185  if (SYNC) ierr = MPI_Barrier(comm);
186  if (TIMING) t1 = MPI_Wtime();
187  ierr = PMPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm);
188  if (TIMING) {TIMER[0] += MPI_Wtime()-t1; TIMER[1] += t1-t0;}
189  return ierr;
190 }
191 
192 #pragma weak MPI_Waitall = mpi_waitall_c
193 #pragma weak MPI_Waitall_ = mpi_waitall_c
194 int mpi_waitall_c(int count, MPI_Request *request, MPI_Status *status)
195 {
196  int ierr;
197  double t0,t1;
198 
199 
200  COUNTER[2]++;
201 
202  if (TIMING) t0 = MPI_Wtime();
203  ierr = PMPI_Waitall(count, request, status);
204  if (TIMING) TIMER[2] += MPI_Wtime()-t0;
205  return ierr;
206 }
207 
208 #pragma weak MPI_Barrier = mpi_barrier_c
209 #pragma weak MPI_Barrier_ = mpi_barrier_c
210 int mpi_barrier_c(MPI_Comm comm)
211 {
212  int ierr;
213  double t0;
214 
215  COUNTER[3]++;
216 
217  if (TIMING) t0 = MPI_Wtime();
218  ierr = PMPI_Barrier(comm);
219  if (TIMING) TIMER[3] += MPI_Wtime()-t0;
220  return ierr;
221 }
222 
223 #pragma weak MPI_Irecv = mpi_irecv_c
224 #pragma weak MPI_Irecv_ = mpi_irecv_c
225 int mpi_irecv_c(void *buf, int count, MPI_Datatype type, int source,
226  int tag, MPI_Comm comm, MPI_Request *request)
227 {
228  int ierr;
229  double t0;
230 
231  COUNTER[4]++;
232 
233  if (TIMING) t0 = MPI_Wtime();
234  ierr = PMPI_Irecv(buf,count,type,source,tag,comm,request);
235  if (TIMING) TIMER[4] += MPI_Wtime()-t0;
236  return ierr;
237 }
238 
239 
240 #pragma weak MPI_Isend = mpi_isend_c
241 #pragma weak MPI_Isend_ = mpi_isend_c
242 int mpi_isend_c(void *buf, int count, MPI_Datatype type, int dest,
243  int tag, MPI_Comm comm, MPI_Request *request)
244 {
245  int ierr;
246  double t0;
247 
248  COUNTER[5]++;
249 
250  if (TIMING) t0 = MPI_Wtime();
251  ierr = PMPI_Isend(buf,count,type,dest,tag,comm,request);
252  if (TIMING) TIMER[5] += MPI_Wtime()-t0;
253  return ierr;
254 }
255 
256 #pragma weak MPI_Recv = mpi_recv_c
257 #pragma weak MPI_Recv_ = mpi_recv_c
258 int mpi_recv_c(void *buf, int count, MPI_Datatype type, int source,
259  int tag, MPI_Comm comm, MPI_Status *status)
260 {
261  int ierr;
262  double t0;
263 
264  COUNTER[6]++;
265 
266  if (TIMING) t0 = MPI_Wtime();
267  ierr = PMPI_Recv(buf,count,type,source,tag,comm,status);
268  if (TIMING) TIMER[6] += MPI_Wtime()-t0;
269  return ierr;
270 }
271 
272 #pragma weak MPI_Send = mpi_send_c
273 #pragma weak MPI_Send_ = mpi_send_c
274 int mpi_send_c(void *buf, int count, MPI_Datatype type, int dest,
275  int tag, MPI_Comm comm)
276 {
277  int ierr;
278  double t0;
279 
280  COUNTER[7]++;
281 
282  if (TIMING) t0 = MPI_Wtime();
283  ierr = PMPI_Send(buf,count,type,dest,tag,comm);
284  if (TIMING) TIMER[7] += MPI_Wtime()-t0;
285  return ierr;
286 }
287 
288 #else
289 
290 void nek_comm_settings(int *sync,int *timing){}
291 void nek_comm_getstat(double *timer, int *counter)
292 {
293  int i;
294  for (i = 0; i < NTIMER; i++) timer[i] = TIMER[i];
295  for (i = 0; i < NCOUNTER; i++) counter[i] = COUNTER[i];
296 }
297 void nek_comm_startstat(void){}
298 
299 #endif
300 
#define nek_comm_settings
Definition: nek_comm.c:34
double TIMER[NTIMER]
Definition: nek_comm.c:41
#define NTIMER
Definition: nek_comm.c:21
void MPI_Comm
Definition: gs_test_old.c:23
int COUNTER[NCOUNTER]
Definition: nek_comm.c:42
int TIMING
Definition: nek_comm.c:40
Definition: comm.h:85
uint count
Definition: xxt.c:179
#define NCOUNTER
Definition: nek_comm.c:22
#define nek_comm_startstat
Definition: nek_comm.c:36
sint datatype
Definition: gs_test_old.c:29
for i
Definition: xxt_test.m:74
int SYNC
Definition: nek_comm.c:39
int MPI_Fint
Definition: comm.h:71
#define nek_comm_getstat
Definition: nek_comm.c:35