xref: /petsc/src/benchmarks/streams/MPIVersion.c (revision 2ff79c18c26c94ed8cb599682f680f231dca6444)
1 /*
2   An adaption of the Stream benchmark for MPI
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 <petscsys.h>
10 
11 #define NTIMESINNER 1
12 //#define N 2*4*20000000
13 #define N 80000000
14 //#define N 1200000
15 //#define N 120000
16 #define NTIMES 50
17 #define OFFSET 0
18 
19 static double a[N+OFFSET],b[N+OFFSET],c[N+OFFSET];
20 static double mintime = FLT_MAX;
21 static double bytes = 3 * sizeof(double) * N;
22 
23 int main(int argc,char **args)
24 {
25   const double scalar = 3.0;
26   double       times[NTIMES],rate;
27   PetscMPIInt  rank,size,resultlen;
28   char         hostname[MPI_MAX_PROCESSOR_NAME] = {0};
29   PetscInt     n = PETSC_DECIDE,NN;
30 
31   PetscCall(PetscInitialize(&argc,&args,NULL,NULL));
32   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD,&rank));
33   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&size));
34 
35   PetscCallMPI(MPI_Get_processor_name(hostname,&resultlen));(void)resultlen;
36   if (rank) PetscCallMPI(MPI_Send(hostname,MPI_MAX_PROCESSOR_NAME,MPI_CHAR,0,0,MPI_COMM_WORLD));
37   else {
38     // printf("Rank %d host %s\n",0,hostname);
39     for (int j = 1; j < size; ++j) {
40       PetscCallMPI(MPI_Recv(hostname,MPI_MAX_PROCESSOR_NAME,MPI_CHAR,j,0,MPI_COMM_WORLD,MPI_STATUS_IGNORE));
41       // printf("Rank %d host %s\n",j,hostname);
42     }
43   }
44 
45   NN = N;
46   PetscCall(PetscSplitOwnership(MPI_COMM_WORLD,&n,&NN));
47   for (int j = 0; j < n; ++j) {
48     a[j] = 1.0;
49     b[j] = 2.0;
50     c[j] = 3.0;
51   }
52 
53   /*   --- MAIN LOOP --- repeat test cases NTIMES times --- */
54   for (int k = 0; k < NTIMES; ++k) {
55     PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD));
56     // Do not include barrier in the timed region
57     times[k] = MPI_Wtime();
58     for (int l=0; l<NTIMESINNER; l++){
59       for (register int j = 0; j < n; j++) a[j] = b[j]+scalar*c[j];
60       if (size == 65) printf("never printed %g\n",a[11]);
61     }
62     //   PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD));
63     times[k] = MPI_Wtime() - times[k];
64   }
65   // use maximum time over all MPI processes
66   PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE,times,NTIMES,MPI_DOUBLE,MPI_MAX,MPI_COMM_WORLD));
67   for (int k = 1; k < NTIMES; ++k) {   /* note -- skip first iteration */
68     mintime = PetscMin(mintime,times[k]);
69   }
70   rate = 1.0E-06 * bytes*NTIMESINNER/mintime;
71 
72   if (rank == 0) {
73     FILE *fd;
74 
75     if (size != 1) {
76       double prate;
77 
78       fd = fopen("flops","r");
79       fscanf(fd,"%lg",&prate);
80       fclose(fd);
81       printf("%d %11.4f   Rate (MB/s) %g \n",size,rate,rate/prate);
82     } else {
83       fd = fopen("flops","w");
84       fprintf(fd,"%g\n",rate);
85       fclose(fd);
86       printf("%d %11.4f   Rate (MB/s) 1\n",size,rate);
87     }
88   }
89   PetscCall(PetscFinalize());
90   return 0;
91 }
92