xref: /petsc/src/benchmarks/streams/BasicVersion.c (revision db4deed73f06ae556a0f2e02a0fd6e016e156849)
15d28107eSBarry Smith 
25d28107eSBarry Smith #include <sys/time.h>
35d28107eSBarry Smith /* int gettimeofday(struct timeval *tp, struct timezone *tzp); */
45d28107eSBarry Smith 
55d28107eSBarry Smith double second()
65d28107eSBarry Smith {
75d28107eSBarry Smith /* struct timeval { long  tv_sec;
85d28107eSBarry Smith                     long  tv_usec; };
95d28107eSBarry Smith 
105d28107eSBarry Smith struct timezone { int tz_minuteswest;
115d28107eSBarry Smith                   int tz_dsttime; }; */
125d28107eSBarry Smith 
135d28107eSBarry Smith   struct timeval tp;
145d28107eSBarry Smith   struct timezone tzp;
155d28107eSBarry Smith   int i;
165d28107eSBarry Smith 
175d28107eSBarry Smith   i = gettimeofday(&tp,&tzp);
185d28107eSBarry Smith   return ( (double) tp.tv_sec + (double) tp.tv_usec * 1.e-6 );
195d28107eSBarry Smith }
205d28107eSBarry Smith # include <stdio.h>
215d28107eSBarry Smith # include <math.h>
225d28107eSBarry Smith # include <limits.h>
230d04baf8SBarry Smith # include <float.h>
245d28107eSBarry Smith # include <sys/time.h>
255d28107eSBarry Smith 
265d28107eSBarry Smith /*
275d28107eSBarry Smith * Program: Stream
285d28107eSBarry Smith * Programmer: Joe R. Zagar
295d28107eSBarry Smith * Revision: 4.0-BETA, October 24, 1995
305d28107eSBarry Smith * Original code developed by John D. McCalpin
315d28107eSBarry Smith *
325d28107eSBarry Smith * This program measures memory transfer rates in MB/s for simple
335d28107eSBarry Smith * computational kernels coded in C.  These numbers reveal the quality
345d28107eSBarry Smith * of code generation for simple uncacheable kernels as well as showing
355d28107eSBarry Smith * the cost of floating-point operations relative to memory accesses.
365d28107eSBarry Smith *
375d28107eSBarry Smith * INSTRUCTIONS:
385d28107eSBarry Smith *
395d28107eSBarry Smith *       1) Stream requires a good bit of memory to run.  Adjust the
405d28107eSBarry Smith *          value of 'N' (below) to give a 'timing calibration' of
415d28107eSBarry Smith *          at least 20 clock-ticks.  This will provide rate estimates
425d28107eSBarry Smith *          that should be good to about 5% precision.
435d28107eSBarry Smith */
445d28107eSBarry Smith 
455d28107eSBarry Smith # define N      2000000
465d28107eSBarry Smith # define NTIMES 50
475d28107eSBarry Smith # define OFFSET 0
485d28107eSBarry Smith 
495d28107eSBarry Smith /*
505d28107eSBarry Smith *      3) Compile the code with full optimization.  Many compilers
515d28107eSBarry Smith *         generate unreasonably bad code before the optimizer tightens
525d28107eSBarry Smith *         things up.  If the results are unreasonably good, on the
535d28107eSBarry Smith *         other hand, the optimizer might be too smart for me!
545d28107eSBarry Smith *
555d28107eSBarry Smith *         Try compiling with:
565d28107eSBarry Smith *               cc -O stream_d.c second.c -o stream_d -lm
575d28107eSBarry Smith *
585d28107eSBarry Smith *         This is known to work on Cray, SGI, IBM, and Sun machines.
595d28107eSBarry Smith *
605d28107eSBarry Smith *
615d28107eSBarry Smith *      4) Mail the results to mccalpin@cs.virginia.edu
625d28107eSBarry Smith *         Be sure to include:
635d28107eSBarry Smith *              a) computer hardware model number and software revision
645d28107eSBarry Smith *              b) the compiler flags
655d28107eSBarry Smith *              c) all of the output from the test case.
665d28107eSBarry Smith * Thanks!
675d28107eSBarry Smith *
685d28107eSBarry Smith */
695d28107eSBarry Smith 
705d28107eSBarry Smith # define HLINE "-------------------------------------------------------------\n"
715d28107eSBarry Smith 
725d28107eSBarry Smith # ifndef MIN
735d28107eSBarry Smith # define MIN(x,y) ((x)<(y)?(x):(y))
745d28107eSBarry Smith # endif
755d28107eSBarry Smith # ifndef MAX
765d28107eSBarry Smith # define MAX(x,y) ((x)>(y)?(x):(y))
775d28107eSBarry Smith # endif
785d28107eSBarry Smith 
795d28107eSBarry Smith static double a[N+OFFSET],
805d28107eSBarry Smith               b[N+OFFSET],
815d28107eSBarry Smith               c[N+OFFSET];
825d28107eSBarry Smith /*double *a,*b,*c;*/
835d28107eSBarry Smith 
84df4a11deSBarry Smith static double mintime[4] = {FLT_MAX,FLT_MAX,FLT_MAX,FLT_MAX};
855d28107eSBarry Smith 
86df4a11deSBarry Smith static const char     *label[4] = {"Copy:      ", "Scale:     ", "Add:       ", "Triad:     "};
875d28107eSBarry Smith 
885d28107eSBarry Smith static double   bytes[4] = {
895d28107eSBarry Smith    2 * sizeof(double) * N,
905d28107eSBarry Smith    2 * sizeof(double) * N,
915d28107eSBarry Smith    3 * sizeof(double) * N,
925d28107eSBarry Smith    3 * sizeof(double) * N
935d28107eSBarry Smith    };
945d28107eSBarry Smith 
955d28107eSBarry Smith extern double second();
965d28107eSBarry Smith 
97c6db04a5SJed Brown #include <mpi.h>
9801a79839SBarry Smith 
9901a79839SBarry Smith int main(int argc,char **args)
1005d28107eSBarry Smith    {
1015d28107eSBarry Smith    int          quantum, checktick();
1025d28107eSBarry Smith    register int j, k;
103df4a11deSBarry Smith    double       scalar, t, times[4][NTIMES],irate[4],rate[4];
10480094aa7SBarry Smith    int          rank,size;
1055d28107eSBarry Smith 
10601a79839SBarry Smith    MPI_Init(&argc,&args);
107df4a11deSBarry Smith    MPI_Comm_rank(MPI_COMM_WORLD,&rank);
10880094aa7SBarry Smith    MPI_Comm_size(MPI_COMM_WORLD,&size);
10980094aa7SBarry Smith    if (!rank)   printf("Number of MPI processes %d\n",size);
110df4a11deSBarry Smith 
1115d28107eSBarry Smith    /* --- SETUP --- determine precision and check timing --- */
1125d28107eSBarry Smith 
113df4a11deSBarry Smith    if (!rank) {
11480094aa7SBarry Smith      /*printf(HLINE);
1155d28107eSBarry Smith      printf("Array size = %d, Offset = %d\n" , N, OFFSET);
116df4a11deSBarry Smith      printf("Total memory required = %.1f MB.\n", (3 * N * BytesPerWord) / 1048576.0);
1175d28107eSBarry Smith      printf("Each test is run %d times, but only\n", NTIMES);
1185d28107eSBarry Smith      printf("the *best* time for each is used.\n");
11980094aa7SBarry Smith       printf(HLINE); */
120df4a11deSBarry Smith    }
1215d28107eSBarry Smith 
1225d28107eSBarry Smith    /* Get initial value for system clock. */
1235d28107eSBarry Smith 
1245d28107eSBarry Smith    /*  a = malloc(N*sizeof(double));
1255d28107eSBarry Smith    b = malloc(N*sizeof(double));
1265d28107eSBarry Smith    c = malloc(N*sizeof(double));*/
1275d28107eSBarry Smith    for (j=0; j<N; j++) {
1285d28107eSBarry Smith         a[j] = 1.0;
1295d28107eSBarry Smith         b[j] = 2.0;
1305d28107eSBarry Smith         c[j] = 0.0;
1315d28107eSBarry Smith         }
1325d28107eSBarry Smith 
133df4a11deSBarry Smith    if (!rank) {
134*db4deed7SKarl Rupp      if  ( (quantum = checktick()) >= 1) ;/* printf("Your clock granularity/precision appears to be %d microseconds.\n", quantum); */
135*db4deed7SKarl Rupp      else ;/* printf("Your clock granularity appears to be less than one microsecond.\n");*/
136df4a11deSBarry Smith    }
1375d28107eSBarry Smith 
1385d28107eSBarry Smith    t = second();
1395d28107eSBarry Smith    for (j = 0; j < N; j++)
1405d28107eSBarry Smith         a[j] = 2.0E0 * a[j];
1415d28107eSBarry Smith    t = 1.0E6 * (second() - t);
1425d28107eSBarry Smith 
143df4a11deSBarry Smith    if (!rank) {
14480094aa7SBarry Smith      /*  printf("Each test below will take on the order of %d microseconds.\n", (int) t  );
1455d28107eSBarry Smith      printf("   (= %d clock ticks)\n", (int) (t/quantum) );
1465d28107eSBarry Smith      printf("Increase the size of the arrays if this shows that\n");
1475d28107eSBarry Smith       printf("you are not getting at least 20 clock ticks per test.\n");
14880094aa7SBarry Smith       printf(HLINE);*/
149df4a11deSBarry Smith    }
1505d28107eSBarry Smith 
1515d28107eSBarry Smith 
1525d28107eSBarry Smith    /*   --- MAIN LOOP --- repeat test cases NTIMES times --- */
1535d28107eSBarry Smith 
1545d28107eSBarry Smith    scalar = 3.0;
1555d28107eSBarry Smith    for (k=0; k<NTIMES; k++)
1565d28107eSBarry Smith         {
15701a79839SBarry Smith    MPI_Barrier(MPI_COMM_WORLD);
1585d28107eSBarry Smith         times[0][k] = second();
159df4a11deSBarry Smith    /* should all these barriers be pulled outside of the time call? */
160106b1c52SJed Brown    MPI_Barrier(MPI_COMM_WORLD);
1615d28107eSBarry Smith         for (j=0; j<N; j++)
1625d28107eSBarry Smith             c[j] = a[j];
163106b1c52SJed Brown    MPI_Barrier(MPI_COMM_WORLD);
1645d28107eSBarry Smith         times[0][k] = second() - times[0][k];
1655d28107eSBarry Smith 
1665d28107eSBarry Smith         times[1][k] = second();
167106b1c52SJed Brown    MPI_Barrier(MPI_COMM_WORLD);
1685d28107eSBarry Smith         for (j=0; j<N; j++)
1695d28107eSBarry Smith             b[j] = scalar*c[j];
170106b1c52SJed Brown    MPI_Barrier(MPI_COMM_WORLD);
1715d28107eSBarry Smith         times[1][k] = second() - times[1][k];
1725d28107eSBarry Smith 
1735d28107eSBarry Smith         times[2][k] = second();
174106b1c52SJed Brown    MPI_Barrier(MPI_COMM_WORLD);
1755d28107eSBarry Smith         for (j=0; j<N; j++)
1765d28107eSBarry Smith             c[j] = a[j]+b[j];
177106b1c52SJed Brown    MPI_Barrier(MPI_COMM_WORLD);
1785d28107eSBarry Smith         times[2][k] = second() - times[2][k];
1795d28107eSBarry Smith 
1805d28107eSBarry Smith         times[3][k] = second();
181106b1c52SJed Brown    MPI_Barrier(MPI_COMM_WORLD);
1825d28107eSBarry Smith         for (j=0; j<N; j++)
1835d28107eSBarry Smith             a[j] = b[j]+scalar*c[j];
184106b1c52SJed Brown    MPI_Barrier(MPI_COMM_WORLD);
1855d28107eSBarry Smith         times[3][k] = second() - times[3][k];
1865d28107eSBarry Smith      }
1875d28107eSBarry Smith 
1885d28107eSBarry Smith    /*   --- SUMMARY --- */
1895d28107eSBarry Smith 
190df4a11deSBarry Smith    for (k=0; k<NTIMES; k++) {
1915d28107eSBarry Smith         for (j=0; j<4; j++) {
192df4a11deSBarry Smith            mintime[j] = MIN(mintime[j], times[j][k]);
193df4a11deSBarry Smith         }
194df4a11deSBarry Smith       }
1955d28107eSBarry Smith 
196df4a11deSBarry Smith    for (j=0; j<4; j++) {
197df4a11deSBarry Smith      irate[j] = 1.0E-06 * bytes[j]/mintime[j];
198df4a11deSBarry Smith    }
199df4a11deSBarry Smith    MPI_Reduce(irate,rate,4,MPI_DOUBLE,MPI_SUM,0,MPI_COMM_WORLD);
200df4a11deSBarry Smith 
201df4a11deSBarry Smith    if (!rank) {
202df4a11deSBarry Smith      printf("Function      Rate (MB/s) \n");
203df4a11deSBarry Smith      for (j=0; j<4; j++) {
204df4a11deSBarry Smith         printf("%s%11.4f\n", label[j],rate[j]);
205df4a11deSBarry Smith      }
2065d28107eSBarry Smith    }
20701a79839SBarry Smith    MPI_Finalize();
2085d28107eSBarry Smith    return 0;
2095d28107eSBarry Smith }
2105d28107eSBarry Smith 
2115d28107eSBarry Smith # define        M        20
2125d28107eSBarry Smith 
2135d28107eSBarry Smith int
2145d28107eSBarry Smith checktick()
2155d28107eSBarry Smith    {
2165d28107eSBarry Smith    int           i, minDelta, Delta;
2175d28107eSBarry Smith    double        t1, t2, timesfound[M];
2185d28107eSBarry Smith 
2195d28107eSBarry Smith /*  Collect a sequence of M unique time values from the system. */
2205d28107eSBarry Smith 
2215d28107eSBarry Smith    for (i = 0; i < M; i++) {
2225d28107eSBarry Smith         t1 = second();
2235d28107eSBarry Smith         while( ((t2=second()) - t1) < 1.0E-6 )
2245d28107eSBarry Smith             ;
2255d28107eSBarry Smith         timesfound[i] = t1 = t2;
2265d28107eSBarry Smith         }
2275d28107eSBarry Smith 
2285d28107eSBarry Smith /*
2295d28107eSBarry Smith * Determine the minimum difference between these M values.
2305d28107eSBarry Smith * This result will be our estimate (in microseconds) for the
2315d28107eSBarry Smith * clock granularity.
2325d28107eSBarry Smith */
2335d28107eSBarry Smith 
2345d28107eSBarry Smith    minDelta = 1000000;
2355d28107eSBarry Smith    for (i = 1; i < M; i++) {
2365d28107eSBarry Smith         Delta = (int)( 1.0E6 * (timesfound[i]-timesfound[i-1]));
2375d28107eSBarry Smith         minDelta = MIN(minDelta, MAX(Delta,0));
2385d28107eSBarry Smith         }
2395d28107eSBarry Smith 
2405d28107eSBarry Smith    return(minDelta);
2415d28107eSBarry Smith    }
2425d28107eSBarry Smith 
243