CoMD
A Mini-app for Co-Design of Classical Molecular Dynamics.
 All Data Structures Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
parallel.c
Go to the documentation of this file.
1 /// \file
2 /// Wrappers for MPI functions. This should be the only compilation
3 /// unit in the code that directly calls MPI functions. To build a pure
4 /// serial version of the code with no MPI, do not define DO_MPI. If
5 /// DO_MPI is not defined then all MPI functionality is replaced with
6 /// equivalent single task behavior.
7 
8 #include "parallel.h"
9 
10 #ifdef DO_MPI
11 #include <mpi.h>
12 #endif
13 
14 #include <stdio.h>
15 #include <time.h>
16 #include <string.h>
17 #include <assert.h>
18 
19 static int myRank = 0;
20 static int nRanks = 1;
21 
22 #ifdef DO_MPI
23 #ifdef SINGLE
24 #define REAL_MPI_TYPE MPI_FLOAT
25 #else
26 #define REAL_MPI_TYPE MPI_DOUBLE
27 #endif
28 
29 #endif
30 
31 int getNRanks()
32 {
33  return nRanks;
34 }
35 
36 int getMyRank()
37 {
38  return myRank;
39 }
40 
41 /// \details
42 /// For now this is just a check for rank 0 but in principle it could be
43 /// more complex. It is also possible to suppress practically all
44 /// output by causing this function to return 0 for all ranks.
45 int printRank()
46 {
47  if (myRank == 0) return 1;
48  return 0;
49 }
50 
51 void timestampBarrier(const char* msg)
52 {
54  if (! printRank())
55  return;
56  time_t t= time(NULL);
57  char* timeString = ctime(&t);
58  timeString[24] = '\0'; // clobber newline
59  fprintf(screenOut, "%s: %s\n", timeString, msg);
60  fflush(screenOut);
61 }
62 
63 void initParallel(int* argc, char*** argv)
64 {
65 #ifdef DO_MPI
66  MPI_Init(argc, argv);
67  MPI_Comm_rank(MPI_COMM_WORLD, &myRank);
68  MPI_Comm_size(MPI_COMM_WORLD, &nRanks);
69 #endif
70 }
71 
73 {
74 #ifdef DO_MPI
75  MPI_Finalize();
76 #endif
77 }
78 
80 {
81 #ifdef DO_MPI
82  MPI_Barrier(MPI_COMM_WORLD);
83 #endif
84 }
85 
86 /// \param [in] sendBuf Data to send.
87 /// \param [in] sendLen Number of bytes to send.
88 /// \param [in] dest Rank in MPI_COMM_WORLD where data will be sent.
89 /// \param [out] recvBuf Received data.
90 /// \param [in] recvLen Maximum number of bytes to receive.
91 /// \param [in] source Rank in MPI_COMM_WORLD from which to receive.
92 /// \return Number of bytes received.
93 int sendReceiveParallel(void* sendBuf, int sendLen, int dest,
94  void* recvBuf, int recvLen, int source)
95 {
96 #ifdef DO_MPI
97  int bytesReceived;
98  MPI_Status status;
99  MPI_Sendrecv(sendBuf, sendLen, MPI_BYTE, dest, 0,
100  recvBuf, recvLen, MPI_BYTE, source, 0,
101  MPI_COMM_WORLD, &status);
102  MPI_Get_count(&status, MPI_BYTE, &bytesReceived);
103 
104  return bytesReceived;
105 #else
106  assert(source == dest);
107  memcpy(recvBuf, sendBuf, sendLen);
108 
109  return sendLen;
110 #endif
111 }
112 
113 void addIntParallel(int* sendBuf, int* recvBuf, int count)
114 {
115 #ifdef DO_MPI
116  MPI_Allreduce(sendBuf, recvBuf, count, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
117 #else
118  for (int ii=0; ii<count; ++ii)
119  recvBuf[ii] = sendBuf[ii];
120 #endif
121 }
122 
123 void addRealParallel(real_t* sendBuf, real_t* recvBuf, int count)
124 {
125 #ifdef DO_MPI
126  MPI_Allreduce(sendBuf, recvBuf, count, REAL_MPI_TYPE, MPI_SUM, MPI_COMM_WORLD);
127 #else
128  for (int ii=0; ii<count; ++ii)
129  recvBuf[ii] = sendBuf[ii];
130 #endif
131 }
132 
133 void addDoubleParallel(double* sendBuf, double* recvBuf, int count)
134 {
135 #ifdef DO_MPI
136  MPI_Allreduce(sendBuf, recvBuf, count, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
137 #else
138  for (int ii=0; ii<count; ++ii)
139  recvBuf[ii] = sendBuf[ii];
140 #endif
141 }
142 
143 void maxIntParallel(int* sendBuf, int* recvBuf, int count)
144 {
145 #ifdef DO_MPI
146  MPI_Allreduce(sendBuf, recvBuf, count, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
147 #else
148  for (int ii=0; ii<count; ++ii)
149  recvBuf[ii] = sendBuf[ii];
150 #endif
151 }
152 
153 
154 void minRankDoubleParallel(RankReduceData* sendBuf, RankReduceData* recvBuf, int count)
155 {
156 #ifdef DO_MPI
157  MPI_Allreduce(sendBuf, recvBuf, count, MPI_DOUBLE_INT, MPI_MINLOC, MPI_COMM_WORLD);
158 #else
159  for (int ii=0; ii<count; ++ii)
160  {
161  recvBuf[ii].val = sendBuf[ii].val;
162  recvBuf[ii].rank = sendBuf[ii].rank;
163  }
164 #endif
165 }
166 
167 void maxRankDoubleParallel(RankReduceData* sendBuf, RankReduceData* recvBuf, int count)
168 {
169 #ifdef DO_MPI
170  MPI_Allreduce(sendBuf, recvBuf, count, MPI_DOUBLE_INT, MPI_MAXLOC, MPI_COMM_WORLD);
171 #else
172  for (int ii=0; ii<count; ++ii)
173  {
174  recvBuf[ii].val = sendBuf[ii].val;
175  recvBuf[ii].rank = sendBuf[ii].rank;
176  }
177 #endif
178 }
179 
180 /// \param [in] count Length of buf in bytes.
181 void bcastParallel(void* buf, int count, int root)
182 {
183 #ifdef DO_MPI
184  MPI_Bcast(buf, count, MPI_BYTE, root, MPI_COMM_WORLD);
185 #endif
186 }
187 
188 int builtWithMpi(void)
189 {
190 #ifdef DO_MPI
191  return 1;
192 #else
193  return 0;
194 #endif
195 }
196 
197