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