Logo AND Algorithmique Numérique Distribuée

Public GIT Repository
some cleaning in the NAS examples
[simgrid.git] / examples / smpi / NAS / EP / ep.c
1 #include <stdlib.h>
2 #include <stdio.h>
3 #include <string.h>
4 #include <math.h>
5
6 #include "mpi.h"
7 #include "npbparams.h"
8
9 #include "simgrid/instr.h" //TRACE_
10
11 #include "randlc.h"
12
13 #ifndef CLASS
14 #define CLASS 'S'
15 #define NUM_PROCS            1                 
16 #endif
17 #define true 1
18 #define false 0
19
20
21 //---NOTE : all the timers function have been modified to
22 //          avoid global timers (privatize these). 
23       // ----------------------- timers ---------------------
24       void timer_clear(double *onetimer) {
25             //elapsed[n] = 0.0;
26             *onetimer = 0.0;
27       }
28
29       void timer_start(double *onetimer) {
30             *onetimer = MPI_Wtime();
31       }
32
33       void timer_stop(int n,double *elapsed,double *start) {
34             double t, now;
35
36             now = MPI_Wtime();
37             t = now - start[n];
38             elapsed[n] += t;
39       }
40
41       double timer_read(int n, double *elapsed) {  /* ok, useless, but jsut to keep function call */
42             return(elapsed[n]);
43       }
44       /********************************************************************
45        *****************            V R A N L C          ******************
46        *****************                                 *****************/           
47       double vranlc(int n, double x, double a, double *y)
48       {
49         int i;
50         long  i246m1=0x00003FFFFFFFFFFF;
51     long  LLx, Lx, La;
52         double d2m46;
53
54 // This doesn't work, because the compiler does the calculation in 32
55 // bits and overflows. No standard way (without f90 stuff) to specify
56 // that the rhs should be done in 64 bit arithmetic.
57 //     parameter(i246m1=2**46-1)
58
59       d2m46=pow(0.5,46);
60
61 // c Note that the v6 compiler on an R8000 does something stupid with
62 // c the above. Using the following instead (or various other things)
63 // c makes the calculation run almost 10 times as fast.
64 //
65 // c     save d2m46
66 // c      data d2m46/0.0d0/
67 // c      if (d2m46 .eq. 0.0d0) then
68 // c         d2m46 = 0.5d0**46
69 // c      endif
70
71       Lx = (long)x;
72       La = (long)a;
73       //fprintf(stdout,("================== Vranlc ================");
74       //fprintf(stdout,("Before Loop: Lx = " + Lx + ", La = " + La);
75   LLx = Lx;
76   for (i=0; i< n; i++) {
77       Lx   = Lx*La & i246m1 ;
78       LLx = Lx;
79       y[i] = d2m46 * (double)LLx;
80       /*
81          if(i == 0) {
82          fprintf(stdout,("After loop 0:");
83          fprintf(stdout,("Lx = " + Lx + ", La = " + La);
84          fprintf(stdout,("d2m46 = " + d2m46);
85          fprintf(stdout,("LLX(Lx) = " + LLX.doubleValue());
86          fprintf(stdout,("Y[0]" + y[0]);
87          }
88        */
89   }
90
91       x = (double)LLx;
92       /*
93       fprintf(stdout,("Change: Lx = " + Lx);
94       fprintf(stdout,("=============End   Vranlc ================");
95       */
96       return x;
97     }
98
99
100
101 //-------------- the core (unique function) -----------
102       void doTest(int argc, char **argv) {
103       double dum[3] = {1.,1.,1.};
104       double x1, x2, sx, sy, tm, an, tt, gc;
105       double Mops;
106       double epsilon=1.0E-8, a = 1220703125., s=271828183.;
107       double t1, t2, t3, t4; 
108       double sx_verify_value, sy_verify_value, sx_err, sy_err;
109
110 #include "npbparams.h"
111       int    mk=16, 
112          // --> set by make : in npbparams.h
113          //m=28, // for CLASS=A
114          //m=30, // for CLASS=B
115          //npm=2, // NPROCS
116          mm = m-mk, 
117          nn = (int)(pow(2,mm)), 
118          nk = (int)(pow(2,mk)), 
119          nq=10, 
120          np, 
121          node, 
122          no_nodes, 
123          i, 
124          ik, 
125          kk, 
126          l, 
127          k, nit, no_large_nodes,
128          np_add, k_offset, j;
129       int    me, nprocs, root=0, dp_type;
130       int verified, 
131           timers_enabled=true;
132       char  size[500]; // mind the size of the string to represent a big number
133
134       //Use in randlc..
135       int KS = 0;
136       double R23, R46, T23, T46;
137
138       double *qq = (double *) malloc (10000*sizeof(double));
139       double *start = (double *) malloc (64*sizeof(double));
140       double *elapsed = (double *) malloc (64*sizeof(double));
141
142       double *x = (double *) malloc (2*nk*sizeof(double));
143       double *q = (double *) malloc (nq*sizeof(double));
144
145       TRACE_smpi_set_category ("start");
146
147       MPI_Init( &argc, &argv );
148       MPI_Comm_size( MPI_COMM_WORLD, &no_nodes);
149       MPI_Comm_rank( MPI_COMM_WORLD, &node);
150
151 #ifdef USE_MPE
152     MPE_Init_log();
153 #endif
154       root = 0;
155       if (node == root ) {
156
157           /*   Because the size of the problem is too large to store in a 32-bit
158            *   integer for some classes, we put it into a string (for printing).
159            *   Have to strip off the decimal point put in there by the floating
160            *   point print statement (internal file)
161            */
162           fprintf(stdout," NAS Parallel Benchmarks 3.2 -- EP Benchmark");
163           sprintf(size,"%d",pow(2,m+1));
164           //size = size.replace('.', ' ');
165           fprintf(stdout," Number of random numbers generated: %s\n",size);
166           fprintf(stdout," Number of active processes: %d\n",no_nodes);
167
168       }
169       verified = false;
170
171       /* c   Compute the number of "batches" of random number pairs generated 
172          c   per processor. Adjust if the number of processors does not evenly 
173          c   divide the total number
174 */
175
176        np = nn / no_nodes;
177        no_large_nodes = nn % no_nodes;
178        if (node < no_large_nodes) np_add = 1;
179        else np_add = 0;
180        np = np + np_add;
181
182        if (np == 0) {
183              fprintf(stdout,"Too many nodes: %d  %d",no_nodes,nn);
184              MPI_Abort(MPI_COMM_WORLD,1);
185              exit(0); 
186        } 
187
188 /* c   Call the random number generator functions and initialize
189    c   the x-array to reduce the effects of paging on the timings.
190    c   Also, call all mathematical functions that are used. Make
191    c   sure these initializations cannot be eliminated as dead code.
192 */
193
194    //call vranlc(0, dum[1], dum[2], dum[3]);
195    // Array indexes start at 1 in Fortran, 0 in Java
196    vranlc(0, dum[0], dum[1], &(dum[2])); 
197
198    dum[0] = randlc(&(dum[1]),&(dum[2]));
199    /////////////////////////////////
200    for (i=0;i<2*nk;i++) {
201        x[i] = -1e99;
202    }
203    Mops = log(sqrt(abs(1))); 
204
205    /*
206       c---------------------------------------------------------------------
207       c    Synchronize before placing time stamp
208       c---------------------------------------------------------------------
209     */
210         MPI_Barrier( MPI_COMM_WORLD );
211
212         TRACE_smpi_set_category ("ep");
213
214         timer_clear(&(elapsed[1]));
215         timer_clear(&(elapsed[2]));
216         timer_clear(&(elapsed[3]));
217         timer_start(&(start[1]));
218         
219         t1 = a;
220   //fprintf(stdout,("(ep.f:160) t1 = " + t1);
221         t1 = vranlc(0, t1, a, x);
222   //fprintf(stdout,("(ep.f:161) t1 = " + t1);
223   
224         
225 /* c   Compute AN = A ^ (2 * NK) (mod 2^46). */
226         
227         t1 = a;
228   //fprintf(stdout,("(ep.f:165) t1 = " + t1);
229         for (i=1; i <= mk+1; i++) {
230                t2 = randlc(&t1, &t1);
231          //fprintf(stdout,("(ep.f:168)[loop i=" + i +"] t1 = " + t1);
232         } 
233         an = t1;
234   //fprintf(stdout,("(ep.f:172) s = " + s);
235         tt = s;
236         gc = 0.;
237         sx = 0.;
238         sy = 0.;
239         for (i=0; i < nq ; i++) {
240                q[i] = 0.;
241         }
242
243 /*
244     Each instance of this loop may be performed independently. We compute
245     the k offsets separately to take into account the fact that some nodes
246     have more numbers to generate than others
247 */
248
249       if (np_add == 1)
250          k_offset = node * np -1;
251       else
252          k_offset = no_large_nodes*(np+1) + (node-no_large_nodes)*np -1;
253      
254       int stop = false;
255       for(k = 1; k <= np; k++) {
256          stop = false;
257          kk = k_offset + k ;
258          t1 = s;
259          //fprintf(stdout,("(ep.f:193) t1 = " + t1);
260          t2 = an;
261
262 //       Find starting seed t1 for this kk.
263
264          for (i=1;i<=100 && !stop;i++) {
265             ik = kk / 2;
266       //fprintf(stdout,("(ep.f:199) ik = " +ik+", kk = " + kk);
267             if (2 * ik != kk)  {
268                 t3 = randlc(&t1, &t2);
269                 //fprintf(stdout,("(ep.f:200) t1= " +t1 );
270             }
271             if (ik==0)
272                 stop = true;
273             else {
274                t3 = randlc(&t2, &t2);
275                kk = ik;
276            }
277          }
278 //       Compute uniform pseudorandom numbers.
279
280          //if (timers_enabled)  timer_start(3);
281    timer_start(&(start[3]));
282          //call vranlc(2 * nk, t1, a, x)  --> t1 and y are modified
283
284   //fprintf(stdout,">>>>>>>>>>>Before vranlc(l.210)<<<<<<<<<<<<<");
285   //fprintf(stdout,"2*nk = " + (2*nk));
286   //fprintf(stdout,"t1 = " + t1);
287   //fprintf(stdout,"a  = " + a);
288   //fprintf(stdout,"x[0] = " + x[0]);
289   //fprintf(stdout,">>>>>>>>>>>>>>>>>>>>>><<<<<<<<<<<<<<<<<<<<<<");
290         
291   t1 = vranlc(2 * nk, t1, a, x);
292
293   //fprintf(stdout,(">>>>>>>>>>>After  Enter vranlc (l.210)<<<<<<");
294   //fprintf(stdout,("2*nk = " + (2*nk));
295   //fprintf(stdout,("t1 = " + t1);
296   //fprintf(stdout,("a  = " + a);
297   //fprintf(stdout,("x[0] = " + x[0]);
298   //fprintf(stdout,(">>>>>>>>>>>>>>>>>>>>>><<<<<<<<<<<<<<<<<<<<<<");
299         
300          //if (timers_enabled)  timer_stop(3);
301    timer_stop(3,elapsed,start);
302
303 /*       Compute Gaussian deviates by acceptance-rejection method and 
304  *       tally counts in concentric square annuli.  This loop is not 
305  *       vectorizable. 
306  */
307          //if (timers_enabled) timer_start(2);
308     timer_start(&(start[2]));
309          for(i=1; i<=nk;i++) {
310             x1 = 2. * x[2*i-2] -1.0;
311             x2 = 2. * x[2*i-1] - 1.0;
312             t1 = x1*x1 + x2*x2;
313             if (t1 <= 1.) {
314                t2   = sqrt(-2. * log(t1) / t1);
315                t3   = (x1 * t2);
316                t4   = (x2 * t2);
317                l    = (int)(abs(t3) > abs(t4) ? abs(t3) : abs(t4));
318                q[l] = q[l] + 1.;
319                sx   = sx + t3;
320                sy   = sy + t4;
321              }
322     /*
323        if(i == 1) {
324                 fprintf(stdout,"x1 = " + x1);
325                 fprintf(stdout,"x2 = " + x2);
326                 fprintf(stdout,"t1 = " + t1);
327                 fprintf(stdout,"t2 = " + t2);
328                 fprintf(stdout,"t3 = " + t3);
329                 fprintf(stdout,"t4 = " + t4);
330                 fprintf(stdout,"l = " + l);
331                 fprintf(stdout,"q[l] = " + q[l]);
332                 fprintf(stdout,"sx = " + sx);
333                 fprintf(stdout,"sy = " + sy);
334        }
335     */
336            }
337          //if (timers_enabled)  timer_stop(2);
338     timer_stop(2,elapsed,start);
339       }
340
341     TRACE_smpi_set_category ("finalize");
342
343       //int MPI_Allreduce(void *sbuf, void *rbuf, int count, MPI_Datatype dtype, MPI_Op op, MPI_Comm comm)   
344   MPI_Allreduce(&sx, x, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
345   sx = x[0]; //FIXME :  x[0] or x[1] => x[0] because fortran starts with 1
346       MPI_Allreduce(&sy, x, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
347       sy = x[0];
348       MPI_Allreduce(q, x, nq, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
349
350       for(i = 0; i < nq; i++) {
351     q[i] = x[i];
352   }
353   for(i = 0; i < nq; i++) {
354     gc += q[i];
355   }
356
357   timer_stop(1,elapsed,start);
358       tm = timer_read(1,elapsed);
359   MPI_Allreduce(&tm, x, 1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
360   tm = x[0];
361
362   if(node == root) {
363     nit = 0;
364     verified = true;
365
366     if(m == 24) {
367              sx_verify_value = -3.247834652034740E3;
368                   sy_verify_value = -6.958407078382297E3;
369     } else if(m == 25) {
370                 sx_verify_value = -2.863319731645753E3;
371       sy_verify_value = -6.320053679109499E3;
372     } else if(m == 28) {
373             sx_verify_value = -4.295875165629892E3;
374       sy_verify_value = -1.580732573678431E4;
375     } else if(m == 30) {
376             sx_verify_value =  4.033815542441498E4;
377                   sy_verify_value = -2.660669192809235E4;
378     } else if(m == 32) {
379                   sx_verify_value =  4.764367927995374E4;
380                      sy_verify_value = -8.084072988043731E4;
381     } else if(m == 36) {
382             sx_verify_value =  1.982481200946593E5;
383             sy_verify_value = -1.020596636361769E5;
384     } else {
385       verified = false;
386     }
387
388     /*
389     fprintf(stdout,("sx        = " + sx);
390     fprintf(stdout,("sx_verify = " + sx_verify_value);
391     fprintf(stdout,("sy        = " + sy);
392     fprintf(stdout,("sy_verify = " + sy_verify_value);
393     */
394     if(verified) {
395       sx_err = abs((sx - sx_verify_value)/sx_verify_value);
396       sy_err = abs((sy - sy_verify_value)/sy_verify_value);
397       /*
398       fprintf(stdout,("sx_err = " + sx_err);
399       fprintf(stdout,("sy_err = " + sx_err);
400       fprintf(stdout,("epsilon= " + epsilon);
401       */
402       verified = ((sx_err < epsilon) && (sy_err < epsilon));
403     }
404
405     Mops = (pow(2.0, m+1))/tm/1000;
406
407     fprintf(stdout,"EP Benchmark Results:\n");
408     fprintf(stdout,"CPU Time=%d\n",tm);
409     fprintf(stdout,"N = 2^%d\n",m);
410     fprintf(stdout,"No. Gaussain Pairs =%d\n",gc);
411     fprintf(stdout,"Sum = %f %ld\n",sx,sy);
412     fprintf(stdout,"Count:");
413     for(i = 0; i < nq; i++) {
414       fprintf(stdout,"%d\t %ld\n",i,q[i]);
415     }
416
417     /*
418     print_results("EP", _class, m+1, 0, 0, nit, npm, no_nodes, tm, Mops,
419         "Random numbers generated", verified, npbversion,
420         compiletime, cs1, cs2, cs3, cs4, cs5, cs6, cs7) */
421     fprintf(stdout,"\nEP Benchmark Completed\n");
422             fprintf(stdout,"Class           = %s\n", _class);
423     fprintf(stdout,"Size            = %s\n", size);
424     fprintf(stdout,"Iteration       = %d\n", nit);
425     fprintf(stdout,"Time in seconds = %f\n",(tm/1000));
426     fprintf(stdout,"Total processes = %d\n",no_nodes);
427     fprintf(stdout,"Mops/s total    = %f\n",Mops);
428     fprintf(stdout,"Mops/s/process  = %f\n", Mops/no_nodes);
429     fprintf(stdout,"Operation type  = Random number generated\n");
430     if(verified) {
431       fprintf(stdout,"Verification    = SUCCESSFUL\n");
432     } else {
433       fprintf(stdout,"Verification    = UNSUCCESSFUL\n");
434     }
435              fprintf(stdout,"Total time:     %f\n",(timer_read(1,elapsed)/1000));
436              fprintf(stdout,"Gaussian pairs: %f\n",(timer_read(2,elapsed)/1000));
437              fprintf(stdout,"Random numbers: %f\n",(timer_read(3,elapsed)/1000));
438          }
439 #ifdef USE_MPE
440     MPE_Finish_log(argv[0]);
441 #endif
442  
443        MPI_Finalize();
444       }
445
446     int main(int argc, char **argv) {
447        doTest(argc,argv);
448     }