19f0612e4SBarry Smith /* 29f0612e4SBarry Smith An adaption of the Stream benchmark for MPI 39f0612e4SBarry Smith Original code developed by John D. McCalpin 49f0612e4SBarry Smith */ 5d3ae85c4SBarry Smith #include <stdio.h> 6d3ae85c4SBarry Smith #include <math.h> 7d3ae85c4SBarry Smith #include <limits.h> 8d3ae85c4SBarry Smith #include <float.h> 95e71baefSBarry Smith #include <petscsys.h> 10d3ae85c4SBarry Smith 119f0612e4SBarry Smith #define NTIMESINNER 1 129f0612e4SBarry Smith //#define N 2*4*20000000 139f0612e4SBarry Smith #define N 80000000 149f0612e4SBarry Smith //#define N 1200000 159f0612e4SBarry Smith //#define N 120000 16d3ae85c4SBarry Smith #define NTIMES 50 17d3ae85c4SBarry Smith #define OFFSET 0 18d3ae85c4SBarry Smith 1951096aa5SJacob Faibussowitsch static double a[N+OFFSET],b[N+OFFSET],c[N+OFFSET]; 209f0612e4SBarry Smith static double mintime = FLT_MAX; 219f0612e4SBarry Smith static double bytes = 3 * sizeof(double) * N; 22d3ae85c4SBarry Smith 23d3ae85c4SBarry Smith int main(int argc,char **args) 24d3ae85c4SBarry Smith { 2551096aa5SJacob Faibussowitsch const double scalar = 3.0; 269f0612e4SBarry Smith double times[NTIMES],rate; 2751096aa5SJacob Faibussowitsch PetscMPIInt rank,size,resultlen; 2851096aa5SJacob Faibussowitsch char hostname[MPI_MAX_PROCESSOR_NAME] = {0}; 299f0612e4SBarry Smith PetscInt n = PETSC_DECIDE,NN; 30d3ae85c4SBarry Smith 31b8abcfdeSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,NULL,NULL)); 3251096aa5SJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD,&rank)); 3351096aa5SJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&size)); 34d3ae85c4SBarry Smith 3551096aa5SJacob Faibussowitsch PetscCallMPI(MPI_Get_processor_name(hostname,&resultlen));(void)resultlen; 3651096aa5SJacob Faibussowitsch if (rank) PetscCallMPI(MPI_Send(hostname,MPI_MAX_PROCESSOR_NAME,MPI_CHAR,0,0,MPI_COMM_WORLD)); 3751096aa5SJacob Faibussowitsch else { 389f0612e4SBarry Smith // printf("Rank %d host %s\n",0,hostname); 3951096aa5SJacob Faibussowitsch for (int j = 1; j < size; ++j) { 4051096aa5SJacob Faibussowitsch PetscCallMPI(MPI_Recv(hostname,MPI_MAX_PROCESSOR_NAME,MPI_CHAR,j,0,MPI_COMM_WORLD,MPI_STATUS_IGNORE)); 419f0612e4SBarry Smith // printf("Rank %d host %s\n",j,hostname); 426b58a888SBarry Smith } 431df1832dSBarry Smith } 44d3ae85c4SBarry Smith 459f0612e4SBarry Smith NN = N; 469f0612e4SBarry Smith PetscCall(PetscSplitOwnership(MPI_COMM_WORLD,&n,&NN)); 479f0612e4SBarry Smith for (int j = 0; j < n; ++j) { 48d3ae85c4SBarry Smith a[j] = 1.0; 49d3ae85c4SBarry Smith b[j] = 2.0; 509f0612e4SBarry Smith c[j] = 3.0; 51d3ae85c4SBarry Smith } 52d3ae85c4SBarry Smith 53d3ae85c4SBarry Smith /* --- MAIN LOOP --- repeat test cases NTIMES times --- */ 5451096aa5SJacob Faibussowitsch for (int k = 0; k < NTIMES; ++k) { 5551096aa5SJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 569f0612e4SBarry Smith // Do not include barrier in the timed region 579f0612e4SBarry Smith times[k] = MPI_Wtime(); 589f0612e4SBarry Smith for (int l=0; l<NTIMESINNER; l++){ 599f0612e4SBarry Smith for (register int j = 0; j < n; j++) a[j] = b[j]+scalar*c[j]; 609f0612e4SBarry Smith if (size == 65) printf("never printed %g\n",a[11]); 61d3ae85c4SBarry Smith } 629f0612e4SBarry Smith // PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 639f0612e4SBarry Smith times[k] = MPI_Wtime() - times[k]; 6451096aa5SJacob Faibussowitsch } 659f0612e4SBarry Smith // use maximum time over all MPI processes 66*458b0db5SMartin Diehl PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE,times,NTIMES,MPI_DOUBLE,MPI_MAX,MPI_COMM_WORLD)); 679f0612e4SBarry Smith for (int k = 1; k < NTIMES; ++k) { /* note -- skip first iteration */ 689f0612e4SBarry Smith mintime = PetscMin(mintime,times[k]); 699f0612e4SBarry Smith } 709f0612e4SBarry Smith rate = 1.0E-06 * bytes*NTIMESINNER/mintime; 71d3ae85c4SBarry Smith 72c5853193SPierre Jolivet if (rank == 0) { 7351096aa5SJacob Faibussowitsch FILE *fd; 7451096aa5SJacob Faibussowitsch 75ed38be93SPierre Jolivet if (size != 1) { 764198fb66SBarry Smith double prate; 7751096aa5SJacob Faibussowitsch 784198fb66SBarry Smith fd = fopen("flops","r"); 794198fb66SBarry Smith fscanf(fd,"%lg",&prate); 804198fb66SBarry Smith fclose(fd); 819f0612e4SBarry Smith printf("%d %11.4f Rate (MB/s) %g \n",size,rate,rate/prate); 8251096aa5SJacob Faibussowitsch } else { 8351096aa5SJacob Faibussowitsch fd = fopen("flops","w"); 849f0612e4SBarry Smith fprintf(fd,"%g\n",rate); 8551096aa5SJacob Faibussowitsch fclose(fd); 869f0612e4SBarry Smith printf("%d %11.4f Rate (MB/s) 1\n",size,rate); 874198fb66SBarry Smith } 88d3ae85c4SBarry Smith } 8951096aa5SJacob Faibussowitsch PetscCall(PetscFinalize()); 90d3ae85c4SBarry Smith return 0; 91d3ae85c4SBarry Smith } 92