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