root / Pi / C / MPI / Pi_MPI.c @ 308
Historique | Voir | Annoter | Télécharger (8,13 ko)
1 | 9 | equemene | //
|
---|---|---|---|
2 | 9 | equemene | // Estimation of Pi using Monte Carlo exploration process
|
3 | 247 | equemene | // Cecill v2 Emmanuel QUEMENER <emmanuel.quemener@gmail.com>
|
4 | 247 | equemene | // gcc -std=c99 -O3 -o Pi_MPI Pi_MPI.c -lm
|
5 | 9 | equemene | |
6 | 56 | equemene | // Needed for gethostname
|
7 | 247 | equemene | #define _DEFAULT_SOURCE
|
8 | 56 | equemene | #include <sys/unistd.h> |
9 | 56 | equemene | |
10 | 9 | equemene | #include <math.h> |
11 | 9 | equemene | #include <stdio.h> |
12 | 9 | equemene | #include <stdlib.h> |
13 | 25 | equemene | #include <limits.h> |
14 | 9 | equemene | #include <mpi.h> |
15 | 126 | equemene | #include <stddef.h> |
16 | 23 | equemene | #include <sys/time.h> |
17 | 23 | equemene | |
18 | 9 | equemene | // Marsaglia RNG very simple implementation
|
19 | 9 | equemene | #define znew ((z=36969*(z&65535)+(z>>16))<<16) |
20 | 9 | equemene | #define wnew ((w=18000*(w&65535)+(w>>16))&65535) |
21 | 9 | equemene | #define MWC (znew+wnew)
|
22 | 9 | equemene | #define SHR3 (jsr=(jsr=(jsr=jsr^(jsr<<17))^(jsr>>13))^(jsr<<5)) |
23 | 9 | equemene | #define CONG (jcong=69069*jcong+1234567) |
24 | 9 | equemene | #define KISS ((MWC^CONG)+SHR3)
|
25 | 9 | equemene | |
26 | 9 | equemene | #define MWCfp MWC * 2.328306435454494e-10f |
27 | 9 | equemene | #define KISSfp KISS * 2.328306435454494e-10f |
28 | 126 | equemene | #define SHR3fp SHR3 * 2.328306435454494e-10f |
29 | 126 | equemene | #define CONGfp CONG * 2.328306435454494e-10f |
30 | 9 | equemene | |
31 | 9 | equemene | #define ITERATIONS 1000000000 |
32 | 9 | equemene | |
33 | 9 | equemene | #ifdef LONG
|
34 | 25 | equemene | #define LENGTH long long |
35 | 9 | equemene | #else
|
36 | 25 | equemene | #define LENGTH int |
37 | 9 | equemene | #endif
|
38 | 9 | equemene | |
39 | 126 | equemene | typedef struct compute_node { |
40 | 126 | equemene | LENGTH inside; |
41 | 126 | equemene | long int useconds; |
42 | 126 | equemene | } result; |
43 | 126 | equemene | |
44 | 56 | equemene | unsigned int rotl(unsigned int value, int shift) { |
45 | 56 | equemene | return (value << shift) | (value >> (sizeof(value) * CHAR_BIT - shift)); |
46 | 56 | equemene | } |
47 | 56 | equemene | |
48 | 56 | equemene | unsigned int rotr(unsigned int value, int shift) { |
49 | 56 | equemene | return (value >> shift) | (value << (sizeof(value) * CHAR_BIT - shift)); |
50 | 56 | equemene | } |
51 | 56 | equemene | |
52 | 9 | equemene | LENGTH MainLoopGlobal(LENGTH iterations,unsigned int seed_w,unsigned int seed_z) |
53 | 9 | equemene | { |
54 | 126 | equemene | #if defined TCONG
|
55 | 126 | equemene | unsigned int jcong=seed_z; |
56 | 126 | equemene | #elif defined TSHR3
|
57 | 126 | equemene | unsigned int jsr=seed_w; |
58 | 126 | equemene | #elif defined TMWC
|
59 | 9 | equemene | unsigned int z=seed_z; |
60 | 9 | equemene | unsigned int w=seed_w; |
61 | 126 | equemene | #elif defined TKISS
|
62 | 126 | equemene | unsigned int jcong=seed_z; |
63 | 126 | equemene | unsigned int jsr=seed_w; |
64 | 126 | equemene | unsigned int z=seed_z; |
65 | 126 | equemene | unsigned int w=seed_w; |
66 | 126 | equemene | #endif
|
67 | 9 | equemene | |
68 | 126 | equemene | LENGTH total=0;
|
69 | 9 | equemene | |
70 | 9 | equemene | for (LENGTH i=0;i<iterations;i++) { |
71 | 9 | equemene | |
72 | 126 | equemene | #if defined TINT32
|
73 | 126 | equemene | #define THEONE 1073741824 |
74 | 126 | equemene | #if defined TCONG
|
75 | 126 | equemene | unsigned int x=CONG>>17 ; |
76 | 126 | equemene | unsigned int y=CONG>>17 ; |
77 | 126 | equemene | #elif defined TSHR3
|
78 | 126 | equemene | unsigned int x=SHR3>>17 ; |
79 | 126 | equemene | unsigned int y=SHR3>>17 ; |
80 | 126 | equemene | #elif defined TMWC
|
81 | 126 | equemene | unsigned int x=MWC>>17 ; |
82 | 126 | equemene | unsigned int y=MWC>>17 ; |
83 | 126 | equemene | #elif defined TKISS
|
84 | 126 | equemene | unsigned int x=KISS>>17 ; |
85 | 126 | equemene | unsigned int y=KISS>>17 ; |
86 | 126 | equemene | #endif
|
87 | 126 | equemene | #elif defined TINT64
|
88 | 126 | equemene | #define THEONE 4611686018427387904 |
89 | 126 | equemene | #if defined TCONG
|
90 | 126 | equemene | unsigned long x=(unsigned long)(CONG>>1) ; |
91 | 126 | equemene | unsigned long y=(unsigned long)(CONG>>1) ; |
92 | 126 | equemene | #elif defined TSHR3
|
93 | 126 | equemene | unsigned long x=(unsigned long)(SHR3>>1) ; |
94 | 126 | equemene | unsigned long y=(unsigned long)(SHR3>>1) ; |
95 | 126 | equemene | #elif defined TMWC
|
96 | 126 | equemene | unsigned long x=(unsigned long)(MWC>>1) ; |
97 | 126 | equemene | unsigned long y=(unsigned long)(MWC>>1) ; |
98 | 126 | equemene | #elif defined TKISS
|
99 | 126 | equemene | unsigned long x=(unsigned long)(KISS>>1) ; |
100 | 126 | equemene | unsigned long y=(unsigned long)(KISS>>1) ; |
101 | 126 | equemene | #endif
|
102 | 126 | equemene | #elif defined TFP32
|
103 | 126 | equemene | #define THEONE 1.0f |
104 | 126 | equemene | #if defined TCONG
|
105 | 126 | equemene | float x=CONGfp ;
|
106 | 126 | equemene | float y=CONGfp ;
|
107 | 126 | equemene | #elif defined TSHR3
|
108 | 126 | equemene | float x=SHR3fp ;
|
109 | 126 | equemene | float y=SHR3fp ;
|
110 | 126 | equemene | #elif defined TMWC
|
111 | 126 | equemene | float x=MWCfp ;
|
112 | 126 | equemene | float y=MWCfp ;
|
113 | 126 | equemene | #elif defined TKISS
|
114 | 126 | equemene | float x=KISSfp ;
|
115 | 126 | equemene | float y=KISSfp ;
|
116 | 126 | equemene | #endif
|
117 | 126 | equemene | #elif defined TFP64
|
118 | 126 | equemene | #define THEONE 1.0f |
119 | 126 | equemene | #if defined TCONG
|
120 | 126 | equemene | double x=(double)CONGfp ; |
121 | 126 | equemene | double y=(double)CONGfp ; |
122 | 126 | equemene | #elif defined TSHR3
|
123 | 126 | equemene | double x=(double)SHR3fp ; |
124 | 126 | equemene | double y=(double)SHR3fp ; |
125 | 126 | equemene | #elif defined TMWC
|
126 | 126 | equemene | double x=(double)MWCfp ; |
127 | 126 | equemene | double y=(double)MWCfp ; |
128 | 126 | equemene | #elif defined TKISS
|
129 | 126 | equemene | double x=(double)KISSfp ; |
130 | 126 | equemene | double y=(double)KISSfp ; |
131 | 126 | equemene | #endif
|
132 | 126 | equemene | #endif
|
133 | 9 | equemene | |
134 | 9 | equemene | // Matching test
|
135 | 126 | equemene | unsigned long inside=((x*x+y*y) < THEONE) ? 1:0; |
136 | 9 | equemene | total+=inside; |
137 | 9 | equemene | } |
138 | 126 | equemene | |
139 | 9 | equemene | return(total);
|
140 | 9 | equemene | } |
141 | 9 | equemene | |
142 | 9 | equemene | int main(int argc, char *argv[]) { |
143 | 9 | equemene | |
144 | 56 | equemene | unsigned int seed_z=362436069,seed_w=52128862; |
145 | 126 | equemene | LENGTH iterations=ITERATIONS,inside[8192],insides,part_inside,part_iterations;
|
146 | 247 | equemene | int NumberProcesses,rank,rc,tag=1,i; |
147 | 9 | equemene | |
148 | 56 | equemene | char hostname[128]; |
149 | 56 | equemene | |
150 | 56 | equemene | gethostname(hostname, sizeof hostname);
|
151 | 56 | equemene | |
152 | 23 | equemene | struct timeval start,end;
|
153 | 23 | equemene | long int useconds; |
154 | 23 | equemene | |
155 | 9 | equemene | MPI_Status Stat; |
156 | 80 | equemene | |
157 | 9 | equemene | rc = MPI_Init(&argc,&argv); |
158 | 9 | equemene | if (rc != MPI_SUCCESS) {
|
159 | 9 | equemene | printf ("Error starting MPI program. Terminating.\n");
|
160 | 9 | equemene | MPI_Abort(MPI_COMM_WORLD, rc); |
161 | 9 | equemene | } |
162 | 9 | equemene | |
163 | 247 | equemene | MPI_Comm_size(MPI_COMM_WORLD,&NumberProcesses); |
164 | 9 | equemene | MPI_Comm_rank(MPI_COMM_WORLD,&rank); |
165 | 9 | equemene | |
166 | 126 | equemene | const int nitems=2; |
167 | 126 | equemene | int blocklengths[2] = {1,1}; |
168 | 126 | equemene | |
169 | 126 | equemene | #ifdef LONG
|
170 | 126 | equemene | MPI_Datatype types[2] = {MPI_LONG, MPI_LONG};
|
171 | 126 | equemene | #else
|
172 | 126 | equemene | MPI_Datatype types[2] = {MPI_INT, MPI_LONG};
|
173 | 126 | equemene | #endif
|
174 | 126 | equemene | |
175 | 126 | equemene | MPI_Datatype mpi_result_type; |
176 | 126 | equemene | MPI_Aint offsets[2];
|
177 | 126 | equemene | |
178 | 126 | equemene | offsets[0] = offsetof(result, inside);
|
179 | 126 | equemene | offsets[1] = offsetof(result, useconds);
|
180 | 126 | equemene | |
181 | 126 | equemene | MPI_Type_create_struct(nitems, blocklengths, offsets, types, &mpi_result_type); |
182 | 126 | equemene | MPI_Type_commit(&mpi_result_type); |
183 | 126 | equemene | |
184 | 9 | equemene | if (rank==0) { |
185 | 247 | equemene | |
186 | 247 | equemene | struct timeval tv1,tv2;
|
187 | 247 | equemene | |
188 | 9 | equemene | if (argc > 1) { |
189 | 25 | equemene | iterations=(LENGTH)atoll(argv[1]);
|
190 | 9 | equemene | } |
191 | 9 | equemene | else {
|
192 | 9 | equemene | printf("\n\tPi : Estimate Pi with Monte Carlo exploration\n\n");
|
193 | 9 | equemene | printf("\t\t#1 : number of iterations (default 1 billion)\n\n");
|
194 | 9 | equemene | } |
195 | 56 | equemene | |
196 | 126 | equemene | printf ("\tSizeof int = %lld bytes.\n", (long long)sizeof(int)); |
197 | 126 | equemene | printf ("\tSizeof long = %lld bytes.\n", (long long)sizeof(long)); |
198 | 126 | equemene | printf ("\tSizeof long long = %lld bytes.\n", (long long)sizeof(long long)); |
199 | 56 | equemene | |
200 | 126 | equemene | printf ("\tMax int = %u\n", INT_MAX);
|
201 | 126 | equemene | printf ("\tMax long = %ld\n", LONG_MAX);
|
202 | 126 | equemene | printf ("\tMax long long = %lld\n\n", LLONG_MAX);
|
203 | 56 | equemene | |
204 | 247 | equemene | part_iterations=((iterations%NumberProcesses) == 0) ? iterations/NumberProcesses:iterations/NumberProcesses+1 ; |
205 | 247 | equemene | |
206 | 247 | equemene | gettimeofday(&tv1, NULL);
|
207 | 56 | equemene | |
208 | 9 | equemene | // Split part of code
|
209 | 247 | equemene | for (i=1;i<NumberProcesses;i++) { |
210 | 56 | equemene | |
211 | 25 | equemene | #ifdef LONG
|
212 | 80 | equemene | rc = MPI_Send(&part_iterations, 1, MPI_LONG_LONG, i, tag,
|
213 | 80 | equemene | MPI_COMM_WORLD); |
214 | 25 | equemene | #else
|
215 | 80 | equemene | rc = MPI_Send(&part_iterations, 1, MPI_INT, i, tag,
|
216 | 80 | equemene | MPI_COMM_WORLD); |
217 | 80 | equemene | #endif
|
218 | 9 | equemene | } |
219 | 9 | equemene | |
220 | 23 | equemene | gettimeofday(&start,(struct timezone *)0); |
221 | 56 | equemene | |
222 | 9 | equemene | insides=MainLoopGlobal(part_iterations,seed_w,seed_z); |
223 | 56 | equemene | |
224 | 23 | equemene | gettimeofday(&end,(struct timezone *)0); |
225 | 23 | equemene | useconds=(end.tv_sec-start.tv_sec)*1000000+end.tv_usec-start.tv_usec;
|
226 | 56 | equemene | |
227 | 126 | equemene | printf("\tOn %s with rank %i, find %lld inside in %lu useconds.\n",
|
228 | 56 | equemene | hostname,rank,(long long)insides,useconds); |
229 | 56 | equemene | |
230 | 9 | equemene | // Join part of code
|
231 | 247 | equemene | for (i=1;i<NumberProcesses;i++) { |
232 | 126 | equemene | |
233 | 126 | equemene | result recv; |
234 | 126 | equemene | |
235 | 126 | equemene | rc = MPI_Recv(&recv, 1, mpi_result_type, i, tag,
|
236 | 80 | equemene | MPI_COMM_WORLD, &Stat); |
237 | 126 | equemene | |
238 | 126 | equemene | inside[i]=recv.inside; |
239 | 126 | equemene | useconds=recv.useconds; |
240 | 126 | equemene | |
241 | 126 | equemene | printf("\tReceive from %i, find %lld inside in %lu useconds\n",i,(long long)inside[i],useconds); |
242 | 126 | equemene | |
243 | 56 | equemene | insides+=inside[i]; |
244 | 56 | equemene | } |
245 | 56 | equemene | |
246 | 247 | equemene | float pi=4.*(float)insides/(float)(part_iterations*NumberProcesses); |
247 | 247 | equemene | |
248 | 247 | equemene | gettimeofday(&tv2, NULL);
|
249 | 56 | equemene | |
250 | 247 | equemene | double elapsed=(double)((tv2.tv_sec-tv1.tv_sec) * 1000000L + |
251 | 247 | equemene | (tv2.tv_usec-tv1.tv_usec))/1000000;
|
252 | 247 | equemene | |
253 | 247 | equemene | double itops=(double)(part_iterations*NumberProcesses)/elapsed; |
254 | 247 | equemene | |
255 | 247 | equemene | printf("\nParallelRate %i\nElapsed Time %.2f\nItops %.0f\nLogItops %.2f\n",NumberProcesses,elapsed,itops,log10(itops));
|
256 | 56 | equemene | |
257 | 247 | equemene | LENGTH total=((iterations%NumberProcesses)==0)?iterations:(iterations/NumberProcesses+1)*NumberProcesses; |
258 | 247 | equemene | |
259 | 247 | equemene | printf("Inside/Total %ld %ld\nPi estimation %f\n\n",(long int)insides,(long int)total,(4.*(float)insides/total)); |
260 | 247 | equemene | |
261 | 9 | equemene | } |
262 | 9 | equemene | else
|
263 | 9 | equemene | { |
264 | 9 | equemene | // Receive information from master
|
265 | 25 | equemene | #ifdef LONG
|
266 | 80 | equemene | rc = MPI_Recv(&part_iterations, 1, MPI_LONG_LONG, 0, tag, |
267 | 80 | equemene | MPI_COMM_WORLD, &Stat); |
268 | 25 | equemene | #else
|
269 | 80 | equemene | rc = MPI_Recv(&part_iterations, 1, MPI_INT, 0, tag, |
270 | 80 | equemene | MPI_COMM_WORLD, &Stat); |
271 | 25 | equemene | #endif
|
272 | 9 | equemene | |
273 | 126 | equemene | /* Not such a good idea to print information...
|
274 | 126 | equemene | printf("\tOn %s with rank %i, receive from master %lld\n",
|
275 | 126 | equemene | hostname,rank,(long long)part_iterations); */
|
276 | 56 | equemene | |
277 | 23 | equemene | gettimeofday(&start,(struct timezone *)0); |
278 | 23 | equemene | |
279 | 56 | equemene | part_inside=MainLoopGlobal(part_iterations,rotr(seed_w,rank),rotl(seed_z,rank)); |
280 | 9 | equemene | |
281 | 56 | equemene | gettimeofday(&end,(struct timezone *)0); |
282 | 56 | equemene | useconds=(end.tv_sec-start.tv_sec)*1000000+end.tv_usec-start.tv_usec;
|
283 | 126 | equemene | |
284 | 126 | equemene | /*
|
285 | 56 | equemene | printf("\tOn %s with rank %i, find %lld inside in %lu useconds.\n",
|
286 | 56 | equemene | hostname,rank,(long long)part_inside,useconds);
|
287 | 126 | equemene | */
|
288 | 23 | equemene | |
289 | 126 | equemene | result send; |
290 | 126 | equemene | send.inside=part_inside; |
291 | 126 | equemene | send.useconds=useconds; |
292 | 126 | equemene | |
293 | 126 | equemene | rc = MPI_Send(&send, 1, mpi_result_type, 0, tag, MPI_COMM_WORLD); |
294 | 126 | equemene | |
295 | 9 | equemene | } |
296 | 126 | equemene | |
297 | 126 | equemene | MPI_Type_free(&mpi_result_type); |
298 | 9 | equemene | MPI_Finalize(); |
299 | 247 | equemene | |
300 | 9 | equemene | } |