xref: /petsc/src/benchmarks/streams/OpenMPVersion.c (revision 1404853cd07923df46691f50dfd5de8bf2e82a6a)
1 /*
2   A simplification of the Stream benchmark for OpenMP
3   Original code developed by John D. McCalpin
4 */
5 #include <stdio.h>
6 #include <math.h>
7 #include <limits.h>
8 #include <float.h>
9 #include <sys/time.h>
10 #include <stdlib.h>
11 #include <petscsys.h>
12 
13 //#define N 2*4*20000000
14 #define N 80000000
15 //#define N 1200000
16 //#define N 120000
17 #define NTIMES       50
18 #define OFFSET       0
19 
20 # if !defined(MIN)
21 # define MIN(x,y) ((x)<(y) ? (x) : (y))
22 # endif
23 # if !defined(MAX)
24 # define MAX(x,y) ((x)>(y) ? (x) : (y))
25 # endif
26 
27 static double a[N+OFFSET],b[N+OFFSET],c[N+OFFSET];
28 static double  mintime = FLT_MAX;
29 static double bytes =  3 * sizeof(double) * N;
30 
31 int main(int argc,char **argv)
32 {
33   MPI_Init(&argc,&argv);
34   const static double scalar = 3.0;
35 #pragma omp threadprivate(scalar)
36   double       times[NTIMES],rate;
37   int          size;
38   char         *env;
39   FILE         *fd;
40 
41   env = getenv("OMP_NUM_THREADS");
42   if (!env) env = (char *) "1";
43   sscanf(env,"%d",&size);
44 
45 #pragma omp parallel for schedule(static)
46   for (int j=0; j<N; j++) {
47     a[j] = 1.0;
48     b[j] = 2.0;
49     c[j] = 3.0;
50   }
51 
52   /*  --- MAIN LOOP --- repeat test cases NTIMES times --- */
53   for (int k=0; k<NTIMES; k++)
54   {
55     times[k] = MPI_Wtime();
56     // https://www.openmp.org/wp-content/uploads/OpenMP-API-Specification-5-2.pdf
57     // #pragma omp parallel for  (same performance as below)
58     // #pragma omp parallel for simd schedule(static)  (same performance as below)
59 #pragma omp parallel for schedule(static)
60     for (register int j=0; j<N; j++) a[j] = b[j]+scalar*c[j];
61     times[k] = MPI_Wtime() - times[k];
62   }
63   for (int k=1; k<NTIMES; k++) {  /* note -- skip first iteration */
64       mintime = MIN(mintime, times[k]);
65   }
66 
67   if (size == 65) printf("Never printed %g\n",a[11]);
68   rate = 1.0E-06 * bytes/mintime;
69 
70   if (size == 1) {
71     printf("%d %11.4f   Rate (MB/s) 1\n",size, rate);
72     fd = fopen("flops","w");
73     fprintf(fd,"%g\n",rate);
74     fclose(fd);
75   } else {
76     double prate;
77     fd = fopen("flops","r");
78     fscanf(fd,"%lg",&prate);
79     fclose(fd);
80     printf("%d %11.4f   Rate (MB/s) %g \n", size, rate,rate/prate);
81   }
82   MPI_Finalize();
83   return 0;
84 }
85