xref: /petsc/src/benchmarks/streams/MPIVersion.c (revision 4198fb6675f53a45d2e0dc864f027d60a0670f7b)
1d3ae85c4SBarry Smith 
2d3ae85c4SBarry Smith # include <stdio.h>
3d3ae85c4SBarry Smith # include <math.h>
4d3ae85c4SBarry Smith # include <limits.h>
5d3ae85c4SBarry Smith # include <float.h>
65e71baefSBarry Smith #include <petscsys.h>
7d3ae85c4SBarry Smith 
8d3ae85c4SBarry Smith /*
9d3ae85c4SBarry Smith * Program: Stream
10d3ae85c4SBarry Smith * Programmer: Joe R. Zagar
11d3ae85c4SBarry Smith * Revision: 4.0-BETA, October 24, 1995
12d3ae85c4SBarry Smith * Original code developed by John D. McCalpin
13d3ae85c4SBarry Smith *
14d3ae85c4SBarry Smith * This program measures memory transfer rates in MB/s for simple
15d3ae85c4SBarry Smith * computational kernels coded in C.  These numbers reveal the quality
16d3ae85c4SBarry Smith * of code generation for simple uncacheable kernels as well as showing
17d3ae85c4SBarry Smith * the cost of floating-point operations relative to memory accesses.
18d3ae85c4SBarry Smith *
19d3ae85c4SBarry Smith * INSTRUCTIONS:
20d3ae85c4SBarry Smith *
21d3ae85c4SBarry Smith *       1) Stream requires a good bit of memory to run.  Adjust the
22d3ae85c4SBarry Smith *          value of 'N' (below) to give a 'timing calibration' of
23d3ae85c4SBarry Smith *          at least 20 clock-ticks.  This will provide rate estimates
24d3ae85c4SBarry Smith *          that should be good to about 5% precision.
25d3ae85c4SBarry Smith */
26d3ae85c4SBarry Smith 
27d3ae85c4SBarry Smith # define N      2000000
28d3ae85c4SBarry Smith # define NTIMES 50
29d3ae85c4SBarry Smith # define OFFSET 0
30d3ae85c4SBarry Smith 
31d3ae85c4SBarry Smith /*
32d3ae85c4SBarry Smith *      3) Compile the code with full optimization.  Many compilers
33d3ae85c4SBarry Smith *         generate unreasonably bad code before the optimizer tightens
34d3ae85c4SBarry Smith *         things up.  If the results are unreasonably good, on the
35d3ae85c4SBarry Smith *         other hand, the optimizer might be too smart for me!
36d3ae85c4SBarry Smith *
37d3ae85c4SBarry Smith *         Try compiling with:
38d3ae85c4SBarry Smith *               cc -O stream_d.c second.c -o stream_d -lm
39d3ae85c4SBarry Smith *
40d3ae85c4SBarry Smith *         This is known to work on Cray, SGI, IBM, and Sun machines.
41d3ae85c4SBarry Smith *
42d3ae85c4SBarry Smith *
43d3ae85c4SBarry Smith *      4) Mail the results to mccalpin@cs.virginia.edu
44d3ae85c4SBarry Smith *         Be sure to include:
45d3ae85c4SBarry Smith *              a) computer hardware model number and software revision
46d3ae85c4SBarry Smith *              b) the compiler flags
47d3ae85c4SBarry Smith *              c) all of the output from the test case.
48d3ae85c4SBarry Smith * Thanks!
49d3ae85c4SBarry Smith *
50d3ae85c4SBarry Smith */
51d3ae85c4SBarry Smith 
52d3ae85c4SBarry Smith # define HLINE "-------------------------------------------------------------\n"
53d3ae85c4SBarry Smith 
54d3ae85c4SBarry Smith # ifndef MIN
55d3ae85c4SBarry Smith # define MIN(x,y) ((x)<(y) ? (x) : (y))
56d3ae85c4SBarry Smith # endif
57d3ae85c4SBarry Smith # ifndef MAX
58d3ae85c4SBarry Smith # define MAX(x,y) ((x)>(y) ? (x) : (y))
59d3ae85c4SBarry Smith # endif
60d3ae85c4SBarry Smith 
61d3ae85c4SBarry Smith static double a[N+OFFSET],
62d3ae85c4SBarry Smith               b[N+OFFSET],
63d3ae85c4SBarry Smith               c[N+OFFSET];
64d3ae85c4SBarry Smith /*double *a,*b,*c;*/
65d3ae85c4SBarry Smith 
66d3ae85c4SBarry Smith static double mintime[4] = {FLT_MAX,FLT_MAX,FLT_MAX,FLT_MAX};
67d3ae85c4SBarry Smith 
68d3ae85c4SBarry Smith static double bytes[4] = {
69d3ae85c4SBarry Smith   2 * sizeof(double) * N,
70d3ae85c4SBarry Smith   2 * sizeof(double) * N,
71d3ae85c4SBarry Smith   3 * sizeof(double) * N,
72d3ae85c4SBarry Smith   3 * sizeof(double) * N
73d3ae85c4SBarry Smith };
74d3ae85c4SBarry Smith 
75d3ae85c4SBarry Smith int main(int argc,char **args)
76d3ae85c4SBarry Smith {
77d1d3a73cSBarry Smith   int            quantum, checktick(void);
78d3ae85c4SBarry Smith   register int   j, k;
79d3ae85c4SBarry Smith   double         scalar, t, times[4][NTIMES],irate[4],rate[4];
80d3ae85c4SBarry Smith   int            rank,size,resultlen;
81d3ae85c4SBarry Smith   char           hostname[MPI_MAX_PROCESSOR_NAME];
821df1832dSBarry Smith   MPI_Status     status;
838a4d7553SBarry Smith   int            ierr;
84*4198fb66SBarry Smith   FILE           *fd;
85d3ae85c4SBarry Smith 
865e71baefSBarry Smith   ierr = PetscInitialize(&argc,&args,NULL,NULL);if (ierr) return ierr;
878a4d7553SBarry Smith   ierr = MPI_Comm_rank(MPI_COMM_WORLD,&rank);if (ierr) return ierr;
888a4d7553SBarry Smith   ierr = MPI_Comm_size(MPI_COMM_WORLD,&size);if (ierr) return ierr;
89d3ae85c4SBarry Smith 
906b58a888SBarry Smith   for (j=0; j<MPI_MAX_PROCESSOR_NAME; j++) {
916b58a888SBarry Smith     hostname[j] = 0;
926b58a888SBarry Smith   }
938a4d7553SBarry Smith   ierr = MPI_Get_processor_name(hostname,&resultlen);if (ierr) return ierr;
941df1832dSBarry Smith   if (!rank) {
951df1832dSBarry Smith     for (j=1; j<size; j++) {
968a4d7553SBarry Smith       ierr = MPI_Recv(hostname,MPI_MAX_PROCESSOR_NAME,MPI_CHAR,j,0,MPI_COMM_WORLD,&status);if (ierr) return ierr;
971df1832dSBarry Smith     }
981df1832dSBarry Smith  } else {
998a4d7553SBarry Smith    ierr = MPI_Send(hostname,MPI_MAX_PROCESSOR_NAME,MPI_CHAR,0,0,MPI_COMM_WORLD);if (ierr) return ierr;
100d3ae85c4SBarry Smith  }
101b7250d5dSSatish Balay  ierr = MPI_Barrier(MPI_COMM_WORLD);
102d3ae85c4SBarry Smith 
103d3ae85c4SBarry Smith   /* --- SETUP --- determine precision and check timing --- */
104d3ae85c4SBarry Smith 
105d3ae85c4SBarry Smith   if (!rank) {
106d3ae85c4SBarry Smith     /*printf(HLINE);
107d3ae85c4SBarry Smith     printf("Array size = %d, Offset = %d\n" , N, OFFSET);
108d3ae85c4SBarry Smith     printf("Total memory required = %.1f MB.\n", (3 * N * BytesPerWord) / 1048576.0);
109d3ae85c4SBarry Smith     printf("Each test is run %d times, but only\n", NTIMES);
110d3ae85c4SBarry Smith     printf("the *best* time for each is used.\n");
111d3ae85c4SBarry Smith     printf(HLINE); */
112d3ae85c4SBarry Smith   }
113d3ae85c4SBarry Smith 
114d3ae85c4SBarry Smith   /* Get initial value for system clock. */
115d3ae85c4SBarry Smith 
116d3ae85c4SBarry Smith   /*  a = malloc(N*sizeof(double));
117d3ae85c4SBarry Smith   b = malloc(N*sizeof(double));
118d3ae85c4SBarry Smith   c = malloc(N*sizeof(double));*/
119d3ae85c4SBarry Smith   for (j=0; j<N; j++) {
120d3ae85c4SBarry Smith     a[j] = 1.0;
121d3ae85c4SBarry Smith     b[j] = 2.0;
122d3ae85c4SBarry Smith     c[j] = 0.0;
123d3ae85c4SBarry Smith   }
124d3ae85c4SBarry Smith 
125d3ae85c4SBarry Smith   if (!rank) {
126d3ae85c4SBarry Smith     if  ((quantum = checktick()) >= 1) ; /* printf("Your clock granularity/precision appears to be %d microseconds.\n", quantum); */
127d3ae85c4SBarry Smith     else ; /* printf("Your clock granularity appears to be less than one microsecond.\n");*/
128d3ae85c4SBarry Smith   }
129d3ae85c4SBarry Smith 
13019623ac0SBarry Smith   t = MPI_Wtime();
131d3ae85c4SBarry Smith   for (j = 0; j < N; j++) a[j] = 2.0E0 * a[j];
13219623ac0SBarry Smith   t = 1.0E6 * (MPI_Wtime() - t);
133d3ae85c4SBarry Smith 
134d3ae85c4SBarry Smith   if (!rank) {
135d3ae85c4SBarry Smith     /*  printf("Each test below will take on the order of %d microseconds.\n", (int) t);
136d3ae85c4SBarry Smith     printf("   (= %d clock ticks)\n", (int) (t/quantum));
137d3ae85c4SBarry Smith     printf("Increase the size of the arrays if this shows that\n");
138d3ae85c4SBarry Smith     printf("you are not getting at least 20 clock ticks per test.\n");
139d3ae85c4SBarry Smith     printf(HLINE);*/
140d3ae85c4SBarry Smith   }
141d3ae85c4SBarry Smith 
142d3ae85c4SBarry Smith 
143d3ae85c4SBarry Smith   /*   --- MAIN LOOP --- repeat test cases NTIMES times --- */
144d3ae85c4SBarry Smith 
145d3ae85c4SBarry Smith   scalar = 3.0;
146d3ae85c4SBarry Smith   for (k=0; k<NTIMES; k++)
147d3ae85c4SBarry Smith   {
148b7250d5dSSatish Balay     ierr = MPI_Barrier(MPI_COMM_WORLD);
14919623ac0SBarry Smith     times[0][k] = MPI_Wtime();
150d3ae85c4SBarry Smith     /* should all these barriers be pulled outside of the time call? */
151b7250d5dSSatish Balay     ierr = MPI_Barrier(MPI_COMM_WORLD);
152d3ae85c4SBarry Smith     for (j=0; j<N; j++) c[j] = a[j];
153b7250d5dSSatish Balay     ierr = MPI_Barrier(MPI_COMM_WORLD);
15419623ac0SBarry Smith     times[0][k] = MPI_Wtime() - times[0][k];
155d3ae85c4SBarry Smith 
15619623ac0SBarry Smith     times[1][k] = MPI_Wtime();
157b7250d5dSSatish Balay     ierr = MPI_Barrier(MPI_COMM_WORLD);
158d3ae85c4SBarry Smith     for (j=0; j<N; j++) b[j] = scalar*c[j];
159b7250d5dSSatish Balay     ierr = MPI_Barrier(MPI_COMM_WORLD);
16019623ac0SBarry Smith     times[1][k] = MPI_Wtime() - times[1][k];
161d3ae85c4SBarry Smith 
16219623ac0SBarry Smith     times[2][k] = MPI_Wtime();
163b7250d5dSSatish Balay     ierr = MPI_Barrier(MPI_COMM_WORLD);
164d3ae85c4SBarry Smith     for (j=0; j<N; j++) c[j] = a[j]+b[j];
165b7250d5dSSatish Balay     ierr = MPI_Barrier(MPI_COMM_WORLD);
16619623ac0SBarry Smith     times[2][k] = MPI_Wtime() - times[2][k];
167d3ae85c4SBarry Smith 
16819623ac0SBarry Smith     times[3][k] = MPI_Wtime();
169b7250d5dSSatish Balay     ierr = MPI_Barrier(MPI_COMM_WORLD);
170d3ae85c4SBarry Smith     for (j=0; j<N; j++) a[j] = b[j]+scalar*c[j];
171b7250d5dSSatish Balay     ierr = MPI_Barrier(MPI_COMM_WORLD);
17219623ac0SBarry Smith     times[3][k] = MPI_Wtime() - times[3][k];
173d3ae85c4SBarry Smith   }
174d3ae85c4SBarry Smith 
175d3ae85c4SBarry Smith   /*   --- SUMMARY --- */
176d3ae85c4SBarry Smith 
177d3ae85c4SBarry Smith   for (k=0; k<NTIMES; k++)
178d3ae85c4SBarry Smith     for (j=0; j<4; j++) mintime[j] = MIN(mintime[j], times[j][k]);
179d3ae85c4SBarry Smith 
180d3ae85c4SBarry Smith   for (j=0; j<4; j++) irate[j] = 1.0E-06 * bytes[j]/mintime[j];
181b7250d5dSSatish Balay   ierr = MPI_Reduce(irate,rate,4,MPI_DOUBLE,MPI_SUM,0,MPI_COMM_WORLD);
182b7250d5dSSatish Balay   if (ierr) printf("Error calling MPI\n");
183d3ae85c4SBarry Smith 
184d3ae85c4SBarry Smith   if (!rank) {
185*4198fb66SBarry Smith     if (size == 1) {
186*4198fb66SBarry Smith       printf("%d %11.4f   Rate (MB/s)\n",size, rate[3]);
187*4198fb66SBarry Smith       fd = fopen("flops","w");
188*4198fb66SBarry Smith       fprintf(fd,"%g\n",rate[3]);
189*4198fb66SBarry Smith       fclose(fd);
190*4198fb66SBarry Smith     } else {
191*4198fb66SBarry Smith       double prate;
192*4198fb66SBarry Smith       fd = fopen("flops","r");
193*4198fb66SBarry Smith       fscanf(fd,"%lg",&prate);
194*4198fb66SBarry Smith       fclose(fd);
195*4198fb66SBarry Smith       printf("%d %11.4f   Rate (MB/s) %g \n", size, rate[3],rate[3]/prate);
196*4198fb66SBarry Smith     }
197d3ae85c4SBarry Smith   }
1985e71baefSBarry Smith   PetscFinalize();
199d3ae85c4SBarry Smith   return 0;
200d3ae85c4SBarry Smith }
201d3ae85c4SBarry Smith 
202d3ae85c4SBarry Smith # define        M        20
203d3ae85c4SBarry Smith 
204d1d3a73cSBarry Smith int checktick(void)
205d3ae85c4SBarry Smith {
206d3ae85c4SBarry Smith   int    i, minDelta, Delta;
207d3ae85c4SBarry Smith   double t1, t2, timesfound[M];
208d3ae85c4SBarry Smith 
209d3ae85c4SBarry Smith /*  Collect a sequence of M unique time values from the system. */
210d3ae85c4SBarry Smith 
211d3ae85c4SBarry Smith   for (i = 0; i < M; i++) {
21219623ac0SBarry Smith     t1 = MPI_Wtime();
21319623ac0SBarry Smith     while (((t2=MPI_Wtime()) - t1) < 1.0E-6) ;
214d3ae85c4SBarry Smith     timesfound[i] = t1 = t2;
215d3ae85c4SBarry Smith   }
216d3ae85c4SBarry Smith 
217d3ae85c4SBarry Smith /*
218d3ae85c4SBarry Smith * Determine the minimum difference between these M values.
219d3ae85c4SBarry Smith * This result will be our estimate (in microseconds) for the
220d3ae85c4SBarry Smith * clock granularity.
221d3ae85c4SBarry Smith */
222d3ae85c4SBarry Smith 
223d3ae85c4SBarry Smith   minDelta = 1000000;
224d3ae85c4SBarry Smith   for (i = 1; i < M; i++) {
225d3ae85c4SBarry Smith     Delta    = (int)(1.0E6 * (timesfound[i]-timesfound[i-1]));
226d3ae85c4SBarry Smith     minDelta = MIN(minDelta, MAX(Delta,0));
227d3ae85c4SBarry Smith   }
228d3ae85c4SBarry Smith 
229d3ae85c4SBarry Smith   return(minDelta);
230d3ae85c4SBarry Smith }
231d3ae85c4SBarry Smith 
232