xref: /petsc/src/benchmarks/streams/MPIVersion.c (revision d3ae85c4231eeb9c96a2cdde9142a652e423115c)
1*d3ae85c4SBarry Smith 
2*d3ae85c4SBarry Smith #include <sys/time.h>
3*d3ae85c4SBarry Smith /* int gettimeofday(struct timeval *tp, struct timezone *tzp); */
4*d3ae85c4SBarry Smith 
5*d3ae85c4SBarry Smith double second()
6*d3ae85c4SBarry Smith {
7*d3ae85c4SBarry Smith /* struct timeval { long  tv_sec;
8*d3ae85c4SBarry Smith                     long  tv_usec; };
9*d3ae85c4SBarry Smith 
10*d3ae85c4SBarry Smith struct timezone { int tz_minuteswest;
11*d3ae85c4SBarry Smith                   int tz_dsttime; }; */
12*d3ae85c4SBarry Smith 
13*d3ae85c4SBarry Smith   struct timeval  tp;
14*d3ae85c4SBarry Smith   struct timezone tzp;
15*d3ae85c4SBarry Smith   int             i;
16*d3ae85c4SBarry Smith 
17*d3ae85c4SBarry Smith   i = gettimeofday(&tp,&tzp);
18*d3ae85c4SBarry Smith   return ((double) tp.tv_sec + (double) tp.tv_usec * 1.e-6);
19*d3ae85c4SBarry Smith }
20*d3ae85c4SBarry Smith # include <stdio.h>
21*d3ae85c4SBarry Smith # include <math.h>
22*d3ae85c4SBarry Smith # include <limits.h>
23*d3ae85c4SBarry Smith # include <float.h>
24*d3ae85c4SBarry Smith # include <sys/time.h>
25*d3ae85c4SBarry Smith 
26*d3ae85c4SBarry Smith /*
27*d3ae85c4SBarry Smith * Program: Stream
28*d3ae85c4SBarry Smith * Programmer: Joe R. Zagar
29*d3ae85c4SBarry Smith * Revision: 4.0-BETA, October 24, 1995
30*d3ae85c4SBarry Smith * Original code developed by John D. McCalpin
31*d3ae85c4SBarry Smith *
32*d3ae85c4SBarry Smith * This program measures memory transfer rates in MB/s for simple
33*d3ae85c4SBarry Smith * computational kernels coded in C.  These numbers reveal the quality
34*d3ae85c4SBarry Smith * of code generation for simple uncacheable kernels as well as showing
35*d3ae85c4SBarry Smith * the cost of floating-point operations relative to memory accesses.
36*d3ae85c4SBarry Smith *
37*d3ae85c4SBarry Smith * INSTRUCTIONS:
38*d3ae85c4SBarry Smith *
39*d3ae85c4SBarry Smith *       1) Stream requires a good bit of memory to run.  Adjust the
40*d3ae85c4SBarry Smith *          value of 'N' (below) to give a 'timing calibration' of
41*d3ae85c4SBarry Smith *          at least 20 clock-ticks.  This will provide rate estimates
42*d3ae85c4SBarry Smith *          that should be good to about 5% precision.
43*d3ae85c4SBarry Smith */
44*d3ae85c4SBarry Smith 
45*d3ae85c4SBarry Smith # define N      2000000
46*d3ae85c4SBarry Smith # define NTIMES 50
47*d3ae85c4SBarry Smith # define OFFSET 0
48*d3ae85c4SBarry Smith 
49*d3ae85c4SBarry Smith /*
50*d3ae85c4SBarry Smith *      3) Compile the code with full optimization.  Many compilers
51*d3ae85c4SBarry Smith *         generate unreasonably bad code before the optimizer tightens
52*d3ae85c4SBarry Smith *         things up.  If the results are unreasonably good, on the
53*d3ae85c4SBarry Smith *         other hand, the optimizer might be too smart for me!
54*d3ae85c4SBarry Smith *
55*d3ae85c4SBarry Smith *         Try compiling with:
56*d3ae85c4SBarry Smith *               cc -O stream_d.c second.c -o stream_d -lm
57*d3ae85c4SBarry Smith *
58*d3ae85c4SBarry Smith *         This is known to work on Cray, SGI, IBM, and Sun machines.
59*d3ae85c4SBarry Smith *
60*d3ae85c4SBarry Smith *
61*d3ae85c4SBarry Smith *      4) Mail the results to mccalpin@cs.virginia.edu
62*d3ae85c4SBarry Smith *         Be sure to include:
63*d3ae85c4SBarry Smith *              a) computer hardware model number and software revision
64*d3ae85c4SBarry Smith *              b) the compiler flags
65*d3ae85c4SBarry Smith *              c) all of the output from the test case.
66*d3ae85c4SBarry Smith * Thanks!
67*d3ae85c4SBarry Smith *
68*d3ae85c4SBarry Smith */
69*d3ae85c4SBarry Smith 
70*d3ae85c4SBarry Smith # define HLINE "-------------------------------------------------------------\n"
71*d3ae85c4SBarry Smith 
72*d3ae85c4SBarry Smith # ifndef MIN
73*d3ae85c4SBarry Smith # define MIN(x,y) ((x)<(y) ? (x) : (y))
74*d3ae85c4SBarry Smith # endif
75*d3ae85c4SBarry Smith # ifndef MAX
76*d3ae85c4SBarry Smith # define MAX(x,y) ((x)>(y) ? (x) : (y))
77*d3ae85c4SBarry Smith # endif
78*d3ae85c4SBarry Smith 
79*d3ae85c4SBarry Smith static double a[N+OFFSET],
80*d3ae85c4SBarry Smith               b[N+OFFSET],
81*d3ae85c4SBarry Smith               c[N+OFFSET];
82*d3ae85c4SBarry Smith /*double *a,*b,*c;*/
83*d3ae85c4SBarry Smith 
84*d3ae85c4SBarry Smith static double mintime[4] = {FLT_MAX,FLT_MAX,FLT_MAX,FLT_MAX};
85*d3ae85c4SBarry Smith 
86*d3ae85c4SBarry Smith static const char *label[4] = {"Copy:      ", "Scale:     ", "Add:       ", "Triad:     "};
87*d3ae85c4SBarry Smith 
88*d3ae85c4SBarry Smith static double bytes[4] = {
89*d3ae85c4SBarry Smith   2 * sizeof(double) * N,
90*d3ae85c4SBarry Smith   2 * sizeof(double) * N,
91*d3ae85c4SBarry Smith   3 * sizeof(double) * N,
92*d3ae85c4SBarry Smith   3 * sizeof(double) * N
93*d3ae85c4SBarry Smith };
94*d3ae85c4SBarry Smith 
95*d3ae85c4SBarry Smith extern double second();
96*d3ae85c4SBarry Smith 
97*d3ae85c4SBarry Smith #include <mpi.h>
98*d3ae85c4SBarry Smith 
99*d3ae85c4SBarry Smith int main(int argc,char **args)
100*d3ae85c4SBarry Smith {
101*d3ae85c4SBarry Smith   int          quantum, checktick();
102*d3ae85c4SBarry Smith   register int j, k;
103*d3ae85c4SBarry Smith   double       scalar, t, times[4][NTIMES],irate[4],rate[4];
104*d3ae85c4SBarry Smith   int          rank,size,resultlen;
105*d3ae85c4SBarry Smith   char         hostname[MPI_MAX_PROCESSOR_NAME];
106*d3ae85c4SBarry Smith 
107*d3ae85c4SBarry Smith   MPI_Init(&argc,&args);
108*d3ae85c4SBarry Smith   MPI_Comm_rank(MPI_COMM_WORLD,&rank);
109*d3ae85c4SBarry Smith   MPI_Comm_size(MPI_COMM_WORLD,&size);
110*d3ae85c4SBarry Smith   if (!rank) printf("Number of MPI processes %d\n",size);
111*d3ae85c4SBarry Smith 
112*d3ae85c4SBarry Smith   for (j=0; j<size; j++) {
113*d3ae85c4SBarry Smith     if (rank == j) {
114*d3ae85c4SBarry Smith       MPI_Get_processor_name(hostname,&resultlen);
115*d3ae85c4SBarry Smith       printf("Process %d %s\n",rank,hostname);
116*d3ae85c4SBarry Smith       fflush(stdout);
117*d3ae85c4SBarry Smith     }
118*d3ae85c4SBarry Smith     MPI_Barrier(MPI_COMM_WORLD);
119*d3ae85c4SBarry Smith   }
120*d3ae85c4SBarry Smith 
121*d3ae85c4SBarry Smith   /* --- SETUP --- determine precision and check timing --- */
122*d3ae85c4SBarry Smith 
123*d3ae85c4SBarry Smith   if (!rank) {
124*d3ae85c4SBarry Smith     /*printf(HLINE);
125*d3ae85c4SBarry Smith     printf("Array size = %d, Offset = %d\n" , N, OFFSET);
126*d3ae85c4SBarry Smith     printf("Total memory required = %.1f MB.\n", (3 * N * BytesPerWord) / 1048576.0);
127*d3ae85c4SBarry Smith     printf("Each test is run %d times, but only\n", NTIMES);
128*d3ae85c4SBarry Smith     printf("the *best* time for each is used.\n");
129*d3ae85c4SBarry Smith     printf(HLINE); */
130*d3ae85c4SBarry Smith   }
131*d3ae85c4SBarry Smith 
132*d3ae85c4SBarry Smith   /* Get initial value for system clock. */
133*d3ae85c4SBarry Smith 
134*d3ae85c4SBarry Smith   /*  a = malloc(N*sizeof(double));
135*d3ae85c4SBarry Smith   b = malloc(N*sizeof(double));
136*d3ae85c4SBarry Smith   c = malloc(N*sizeof(double));*/
137*d3ae85c4SBarry Smith   for (j=0; j<N; j++) {
138*d3ae85c4SBarry Smith     a[j] = 1.0;
139*d3ae85c4SBarry Smith     b[j] = 2.0;
140*d3ae85c4SBarry Smith     c[j] = 0.0;
141*d3ae85c4SBarry Smith   }
142*d3ae85c4SBarry Smith 
143*d3ae85c4SBarry Smith   if (!rank) {
144*d3ae85c4SBarry Smith     if  ((quantum = checktick()) >= 1) ; /* printf("Your clock granularity/precision appears to be %d microseconds.\n", quantum); */
145*d3ae85c4SBarry Smith     else ; /* printf("Your clock granularity appears to be less than one microsecond.\n");*/
146*d3ae85c4SBarry Smith   }
147*d3ae85c4SBarry Smith 
148*d3ae85c4SBarry Smith   t = second();
149*d3ae85c4SBarry Smith   for (j = 0; j < N; j++) a[j] = 2.0E0 * a[j];
150*d3ae85c4SBarry Smith   t = 1.0E6 * (second() - t);
151*d3ae85c4SBarry Smith 
152*d3ae85c4SBarry Smith   if (!rank) {
153*d3ae85c4SBarry Smith     /*  printf("Each test below will take on the order of %d microseconds.\n", (int) t);
154*d3ae85c4SBarry Smith     printf("   (= %d clock ticks)\n", (int) (t/quantum));
155*d3ae85c4SBarry Smith     printf("Increase the size of the arrays if this shows that\n");
156*d3ae85c4SBarry Smith     printf("you are not getting at least 20 clock ticks per test.\n");
157*d3ae85c4SBarry Smith     printf(HLINE);*/
158*d3ae85c4SBarry Smith   }
159*d3ae85c4SBarry Smith 
160*d3ae85c4SBarry Smith 
161*d3ae85c4SBarry Smith   /*   --- MAIN LOOP --- repeat test cases NTIMES times --- */
162*d3ae85c4SBarry Smith 
163*d3ae85c4SBarry Smith   scalar = 3.0;
164*d3ae85c4SBarry Smith   for (k=0; k<NTIMES; k++)
165*d3ae85c4SBarry Smith   {
166*d3ae85c4SBarry Smith     MPI_Barrier(MPI_COMM_WORLD);
167*d3ae85c4SBarry Smith     times[0][k] = second();
168*d3ae85c4SBarry Smith     /* should all these barriers be pulled outside of the time call? */
169*d3ae85c4SBarry Smith     MPI_Barrier(MPI_COMM_WORLD);
170*d3ae85c4SBarry Smith     for (j=0; j<N; j++) c[j] = a[j];
171*d3ae85c4SBarry Smith     MPI_Barrier(MPI_COMM_WORLD);
172*d3ae85c4SBarry Smith     times[0][k] = second() - times[0][k];
173*d3ae85c4SBarry Smith 
174*d3ae85c4SBarry Smith     times[1][k] = second();
175*d3ae85c4SBarry Smith     MPI_Barrier(MPI_COMM_WORLD);
176*d3ae85c4SBarry Smith     for (j=0; j<N; j++) b[j] = scalar*c[j];
177*d3ae85c4SBarry Smith     MPI_Barrier(MPI_COMM_WORLD);
178*d3ae85c4SBarry Smith     times[1][k] = second() - times[1][k];
179*d3ae85c4SBarry Smith 
180*d3ae85c4SBarry Smith     times[2][k] = second();
181*d3ae85c4SBarry Smith     MPI_Barrier(MPI_COMM_WORLD);
182*d3ae85c4SBarry Smith     for (j=0; j<N; j++) c[j] = a[j]+b[j];
183*d3ae85c4SBarry Smith     MPI_Barrier(MPI_COMM_WORLD);
184*d3ae85c4SBarry Smith     times[2][k] = second() - times[2][k];
185*d3ae85c4SBarry Smith 
186*d3ae85c4SBarry Smith     times[3][k] = second();
187*d3ae85c4SBarry Smith     MPI_Barrier(MPI_COMM_WORLD);
188*d3ae85c4SBarry Smith     for (j=0; j<N; j++) a[j] = b[j]+scalar*c[j];
189*d3ae85c4SBarry Smith     MPI_Barrier(MPI_COMM_WORLD);
190*d3ae85c4SBarry Smith     times[3][k] = second() - times[3][k];
191*d3ae85c4SBarry Smith   }
192*d3ae85c4SBarry Smith 
193*d3ae85c4SBarry Smith   /*   --- SUMMARY --- */
194*d3ae85c4SBarry Smith 
195*d3ae85c4SBarry Smith   for (k=0; k<NTIMES; k++)
196*d3ae85c4SBarry Smith     for (j=0; j<4; j++) mintime[j] = MIN(mintime[j], times[j][k]);
197*d3ae85c4SBarry Smith 
198*d3ae85c4SBarry Smith   for (j=0; j<4; j++) irate[j] = 1.0E-06 * bytes[j]/mintime[j];
199*d3ae85c4SBarry Smith   MPI_Reduce(irate,rate,4,MPI_DOUBLE,MPI_SUM,0,MPI_COMM_WORLD);
200*d3ae85c4SBarry Smith 
201*d3ae85c4SBarry Smith   if (!rank) {
202*d3ae85c4SBarry Smith     printf("Function      Rate (MB/s) \n");
203*d3ae85c4SBarry Smith     for (j=0; j<4; j++) printf("%s%11.4f\n", label[j],rate[j]);
204*d3ae85c4SBarry Smith   }
205*d3ae85c4SBarry Smith   MPI_Finalize();
206*d3ae85c4SBarry Smith   return 0;
207*d3ae85c4SBarry Smith }
208*d3ae85c4SBarry Smith 
209*d3ae85c4SBarry Smith # define        M        20
210*d3ae85c4SBarry Smith 
211*d3ae85c4SBarry Smith int checktick()
212*d3ae85c4SBarry Smith {
213*d3ae85c4SBarry Smith   int    i, minDelta, Delta;
214*d3ae85c4SBarry Smith   double t1, t2, timesfound[M];
215*d3ae85c4SBarry Smith 
216*d3ae85c4SBarry Smith /*  Collect a sequence of M unique time values from the system. */
217*d3ae85c4SBarry Smith 
218*d3ae85c4SBarry Smith   for (i = 0; i < M; i++) {
219*d3ae85c4SBarry Smith     t1 = second();
220*d3ae85c4SBarry Smith     while (((t2=second()) - t1) < 1.0E-6) ;
221*d3ae85c4SBarry Smith     timesfound[i] = t1 = t2;
222*d3ae85c4SBarry Smith   }
223*d3ae85c4SBarry Smith 
224*d3ae85c4SBarry Smith /*
225*d3ae85c4SBarry Smith * Determine the minimum difference between these M values.
226*d3ae85c4SBarry Smith * This result will be our estimate (in microseconds) for the
227*d3ae85c4SBarry Smith * clock granularity.
228*d3ae85c4SBarry Smith */
229*d3ae85c4SBarry Smith 
230*d3ae85c4SBarry Smith   minDelta = 1000000;
231*d3ae85c4SBarry Smith   for (i = 1; i < M; i++) {
232*d3ae85c4SBarry Smith     Delta    = (int)(1.0E6 * (timesfound[i]-timesfound[i-1]));
233*d3ae85c4SBarry Smith     minDelta = MIN(minDelta, MAX(Delta,0));
234*d3ae85c4SBarry Smith   }
235*d3ae85c4SBarry Smith 
236*d3ae85c4SBarry Smith   return(minDelta);
237*d3ae85c4SBarry Smith }
238*d3ae85c4SBarry Smith 
239