1403adfb6SMatthew G Knepley /* 2403adfb6SMatthew G Knepley STREAM benchmark implementation in CUDA. 3403adfb6SMatthew G Knepley 4403adfb6SMatthew G Knepley COPY: a(i) = b(i) 5403adfb6SMatthew G Knepley SCALE: a(i) = q*b(i) 6403adfb6SMatthew G Knepley SUM: a(i) = b(i) + c(i) 7403adfb6SMatthew G Knepley TRIAD: a(i) = b(i) + q*c(i) 8403adfb6SMatthew G Knepley 9403adfb6SMatthew G Knepley It measures the memory system on the device. 1019816777SMark The implementation is in double precision with a single option. 11403adfb6SMatthew G Knepley 12403adfb6SMatthew G Knepley Code based on the code developed by John D. McCalpin 13403adfb6SMatthew G Knepley http://www.cs.virginia.edu/stream/FTP/Code/stream.c 14403adfb6SMatthew G Knepley 15403adfb6SMatthew G Knepley Written by: Massimiliano Fatica, NVIDIA Corporation 16403adfb6SMatthew G Knepley Modified by: Douglas Enright (dpephd-nvidia@yahoo.com), 1 December 2010 17403adfb6SMatthew G Knepley Extensive Revisions, 4 December 2010 18403adfb6SMatthew G Knepley Modified for PETSc by: Matthew G. Knepley 14 Aug 2011 19403adfb6SMatthew G Knepley 20403adfb6SMatthew G Knepley User interface motivated by bandwidthTest NVIDIA SDK example. 21403adfb6SMatthew G Knepley */ 2219816777SMark static char help[] = "Double-Precision STREAM Benchmark implementation in CUDA\n Performs Copy, Scale, Add, and Triad double-precision kernels\n\n"; 23403adfb6SMatthew G Knepley 24403adfb6SMatthew G Knepley #include <petscconf.h> 25403adfb6SMatthew G Knepley #include <petscsys.h> 26403adfb6SMatthew G Knepley #include <petsctime.h> 275f80ce2aSJacob Faibussowitsch #include <petscdevice.h> 28403adfb6SMatthew G Knepley 2919816777SMark #define N 10000000 30403adfb6SMatthew G Knepley #define NTIMES 10 31403adfb6SMatthew G Knepley 32403adfb6SMatthew G Knepley # ifndef MIN 33403adfb6SMatthew G Knepley # define MIN(x,y) ((x)<(y) ? (x) : (y)) 34403adfb6SMatthew G Knepley # endif 35403adfb6SMatthew G Knepley # ifndef MAX 36403adfb6SMatthew G Knepley # define MAX(x,y) ((x)>(y) ? (x) : (y)) 37403adfb6SMatthew G Knepley # endif 38403adfb6SMatthew G Knepley 39403adfb6SMatthew G Knepley const float flt_eps = 1.192092896e-07f; 40caccb7e3SMatthew G Knepley const double dbl_eps = 2.2204460492503131e-16; 41403adfb6SMatthew G Knepley 42403adfb6SMatthew G Knepley __global__ void set_array(float *a, float value, size_t len) 43403adfb6SMatthew G Knepley { 44403adfb6SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 45403adfb6SMatthew G Knepley while (idx < len) { 46403adfb6SMatthew G Knepley a[idx] = value; 47403adfb6SMatthew G Knepley idx += blockDim.x * gridDim.x; 48403adfb6SMatthew G Knepley } 49403adfb6SMatthew G Knepley } 50403adfb6SMatthew G Knepley 51caccb7e3SMatthew G Knepley __global__ void set_array_double(double *a, double value, size_t len) 52caccb7e3SMatthew G Knepley { 53caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 54caccb7e3SMatthew G Knepley while (idx < len) { 55caccb7e3SMatthew G Knepley a[idx] = value; 56caccb7e3SMatthew G Knepley idx += blockDim.x * gridDim.x; 57caccb7e3SMatthew G Knepley } 58caccb7e3SMatthew G Knepley } 59caccb7e3SMatthew G Knepley 60403adfb6SMatthew G Knepley __global__ void STREAM_Copy(float *a, float *b, size_t len) 61403adfb6SMatthew G Knepley { 62403adfb6SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 63403adfb6SMatthew G Knepley while (idx < len) { 64403adfb6SMatthew G Knepley b[idx] = a[idx]; 65403adfb6SMatthew G Knepley idx += blockDim.x * gridDim.x; 66403adfb6SMatthew G Knepley } 67403adfb6SMatthew G Knepley } 68403adfb6SMatthew G Knepley 69caccb7e3SMatthew G Knepley __global__ void STREAM_Copy_double(double *a, double *b, size_t len) 70caccb7e3SMatthew G Knepley { 71caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 72caccb7e3SMatthew G Knepley while (idx < len) { 73caccb7e3SMatthew G Knepley b[idx] = a[idx]; 74caccb7e3SMatthew G Knepley idx += blockDim.x * gridDim.x; 75caccb7e3SMatthew G Knepley } 76caccb7e3SMatthew G Knepley } 77caccb7e3SMatthew G Knepley 78403adfb6SMatthew G Knepley __global__ void STREAM_Copy_Optimized(float *a, float *b, size_t len) 79403adfb6SMatthew G Knepley { 80403adfb6SMatthew G Knepley /* 81403adfb6SMatthew G Knepley * Ensure size of thread index space is as large as or greater than 82403adfb6SMatthew G Knepley * vector index space else return. 83403adfb6SMatthew G Knepley */ 84403adfb6SMatthew G Knepley if (blockDim.x * gridDim.x < len) return; 85403adfb6SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 86403adfb6SMatthew G Knepley if (idx < len) b[idx] = a[idx]; 87403adfb6SMatthew G Knepley } 88403adfb6SMatthew G Knepley 89caccb7e3SMatthew G Knepley __global__ void STREAM_Copy_Optimized_double(double *a, double *b, size_t len) 90caccb7e3SMatthew G Knepley { 91caccb7e3SMatthew G Knepley /* 92caccb7e3SMatthew G Knepley * Ensure size of thread index space is as large as or greater than 93caccb7e3SMatthew G Knepley * vector index space else return. 94caccb7e3SMatthew G Knepley */ 95caccb7e3SMatthew G Knepley if (blockDim.x * gridDim.x < len) return; 96caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 97caccb7e3SMatthew G Knepley if (idx < len) b[idx] = a[idx]; 98caccb7e3SMatthew G Knepley } 99caccb7e3SMatthew G Knepley 100403adfb6SMatthew G Knepley __global__ void STREAM_Scale(float *a, float *b, float scale, size_t len) 101403adfb6SMatthew G Knepley { 102403adfb6SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 103403adfb6SMatthew G Knepley while (idx < len) { 104403adfb6SMatthew G Knepley b[idx] = scale* a[idx]; 105403adfb6SMatthew G Knepley idx += blockDim.x * gridDim.x; 106403adfb6SMatthew G Knepley } 107403adfb6SMatthew G Knepley } 108403adfb6SMatthew G Knepley 109caccb7e3SMatthew G Knepley __global__ void STREAM_Scale_double(double *a, double *b, double scale, size_t len) 110caccb7e3SMatthew G Knepley { 111caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 112caccb7e3SMatthew G Knepley while (idx < len) { 113caccb7e3SMatthew G Knepley b[idx] = scale* a[idx]; 114caccb7e3SMatthew G Knepley idx += blockDim.x * gridDim.x; 115caccb7e3SMatthew G Knepley } 116caccb7e3SMatthew G Knepley } 117caccb7e3SMatthew G Knepley 118caccb7e3SMatthew G Knepley __global__ void STREAM_Scale_Optimized(float *a, float *b, float scale, size_t len) 119caccb7e3SMatthew G Knepley { 120caccb7e3SMatthew G Knepley /* 121caccb7e3SMatthew G Knepley * Ensure size of thread index space is as large as or greater than 122caccb7e3SMatthew G Knepley * vector index space else return. 123caccb7e3SMatthew G Knepley */ 124caccb7e3SMatthew G Knepley if (blockDim.x * gridDim.x < len) return; 125caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 126caccb7e3SMatthew G Knepley if (idx < len) b[idx] = scale* a[idx]; 127caccb7e3SMatthew G Knepley } 128caccb7e3SMatthew G Knepley 129caccb7e3SMatthew G Knepley __global__ void STREAM_Scale_Optimized_double(double *a, double *b, double scale, size_t len) 130caccb7e3SMatthew G Knepley { 131caccb7e3SMatthew G Knepley /* 132caccb7e3SMatthew G Knepley * Ensure size of thread index space is as large as or greater than 133caccb7e3SMatthew G Knepley * vector index space else return. 134caccb7e3SMatthew G Knepley */ 135caccb7e3SMatthew G Knepley if (blockDim.x * gridDim.x < len) return; 136caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 137caccb7e3SMatthew G Knepley if (idx < len) b[idx] = scale* a[idx]; 138caccb7e3SMatthew G Knepley } 139caccb7e3SMatthew G Knepley 140403adfb6SMatthew G Knepley __global__ void STREAM_Add(float *a, float *b, float *c, size_t len) 141403adfb6SMatthew G Knepley { 142403adfb6SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 143403adfb6SMatthew G Knepley while (idx < len) { 144403adfb6SMatthew G Knepley c[idx] = a[idx]+b[idx]; 145403adfb6SMatthew G Knepley idx += blockDim.x * gridDim.x; 146403adfb6SMatthew G Knepley } 147403adfb6SMatthew G Knepley } 148403adfb6SMatthew G Knepley 149caccb7e3SMatthew G Knepley __global__ void STREAM_Add_double(double *a, double *b, double *c, size_t len) 150caccb7e3SMatthew G Knepley { 151caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 152caccb7e3SMatthew G Knepley while (idx < len) { 153caccb7e3SMatthew G Knepley c[idx] = a[idx]+b[idx]; 154caccb7e3SMatthew G Knepley idx += blockDim.x * gridDim.x; 155caccb7e3SMatthew G Knepley } 156caccb7e3SMatthew G Knepley } 157caccb7e3SMatthew G Knepley 158caccb7e3SMatthew G Knepley __global__ void STREAM_Add_Optimized(float *a, float *b, float *c, size_t len) 159caccb7e3SMatthew G Knepley { 160caccb7e3SMatthew G Knepley /* 161caccb7e3SMatthew G Knepley * Ensure size of thread index space is as large as or greater than 162caccb7e3SMatthew G Knepley * vector index space else return. 163caccb7e3SMatthew G Knepley */ 164caccb7e3SMatthew G Knepley if (blockDim.x * gridDim.x < len) return; 165caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 166caccb7e3SMatthew G Knepley if (idx < len) c[idx] = a[idx]+b[idx]; 167caccb7e3SMatthew G Knepley } 168caccb7e3SMatthew G Knepley 169caccb7e3SMatthew G Knepley __global__ void STREAM_Add_Optimized_double(double *a, double *b, double *c, size_t len) 170caccb7e3SMatthew G Knepley { 171caccb7e3SMatthew G Knepley /* 172caccb7e3SMatthew G Knepley * Ensure size of thread index space is as large as or greater than 173caccb7e3SMatthew G Knepley * vector index space else return. 174caccb7e3SMatthew G Knepley */ 175caccb7e3SMatthew G Knepley if (blockDim.x * gridDim.x < len) return; 176caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 177caccb7e3SMatthew G Knepley if (idx < len) c[idx] = a[idx]+b[idx]; 178caccb7e3SMatthew G Knepley } 179caccb7e3SMatthew G Knepley 180403adfb6SMatthew G Knepley __global__ void STREAM_Triad(float *a, float *b, float *c, float scalar, size_t len) 181403adfb6SMatthew G Knepley { 182403adfb6SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 183403adfb6SMatthew G Knepley while (idx < len) { 184403adfb6SMatthew G Knepley c[idx] = a[idx]+scalar*b[idx]; 185403adfb6SMatthew G Knepley idx += blockDim.x * gridDim.x; 186403adfb6SMatthew G Knepley } 187403adfb6SMatthew G Knepley } 188403adfb6SMatthew G Knepley 189caccb7e3SMatthew G Knepley __global__ void STREAM_Triad_double(double *a, double *b, double *c, double scalar, size_t len) 190caccb7e3SMatthew G Knepley { 191caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 192caccb7e3SMatthew G Knepley while (idx < len) { 193caccb7e3SMatthew G Knepley c[idx] = a[idx]+scalar*b[idx]; 194caccb7e3SMatthew G Knepley idx += blockDim.x * gridDim.x; 195caccb7e3SMatthew G Knepley } 196caccb7e3SMatthew G Knepley } 197caccb7e3SMatthew G Knepley 198caccb7e3SMatthew G Knepley __global__ void STREAM_Triad_Optimized(float *a, float *b, float *c, float scalar, size_t len) 199caccb7e3SMatthew G Knepley { 200caccb7e3SMatthew G Knepley /* 201caccb7e3SMatthew G Knepley * Ensure size of thread index space is as large as or greater than 202caccb7e3SMatthew G Knepley * vector index space else return. 203caccb7e3SMatthew G Knepley */ 204caccb7e3SMatthew G Knepley if (blockDim.x * gridDim.x < len) return; 205caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 206caccb7e3SMatthew G Knepley if (idx < len) c[idx] = a[idx]+scalar*b[idx]; 207caccb7e3SMatthew G Knepley } 208caccb7e3SMatthew G Knepley 209caccb7e3SMatthew G Knepley __global__ void STREAM_Triad_Optimized_double(double *a, double *b, double *c, double scalar, size_t len) 210caccb7e3SMatthew G Knepley { 211caccb7e3SMatthew G Knepley /* 212caccb7e3SMatthew G Knepley * Ensure size of thread index space is as large as or greater than 213caccb7e3SMatthew G Knepley * vector index space else return. 214caccb7e3SMatthew G Knepley */ 215caccb7e3SMatthew G Knepley if (blockDim.x * gridDim.x < len) return; 216caccb7e3SMatthew G Knepley size_t idx = threadIdx.x + blockIdx.x * blockDim.x; 217caccb7e3SMatthew G Knepley if (idx < len) c[idx] = a[idx]+scalar*b[idx]; 218caccb7e3SMatthew G Knepley } 219caccb7e3SMatthew G Knepley 220403adfb6SMatthew G Knepley /* Host side verification routines */ 221a6dfd86eSKarl Rupp bool STREAM_Copy_verify(float *a, float *b, size_t len) 222a6dfd86eSKarl Rupp { 223403adfb6SMatthew G Knepley size_t idx; 224403adfb6SMatthew G Knepley bool bDifferent = false; 225403adfb6SMatthew G Knepley 226403adfb6SMatthew G Knepley for (idx = 0; idx < len && !bDifferent; idx++) { 227403adfb6SMatthew G Knepley float expectedResult = a[idx]; 228403adfb6SMatthew G Knepley float diffResultExpected = (b[idx] - expectedResult); 229403adfb6SMatthew G Knepley float relErrorULPS = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps; 230403adfb6SMatthew G Knepley /* element-wise relative error determination */ 231403adfb6SMatthew G Knepley bDifferent = (relErrorULPS > 2.f); 232403adfb6SMatthew G Knepley } 233403adfb6SMatthew G Knepley 234403adfb6SMatthew G Knepley return bDifferent; 235403adfb6SMatthew G Knepley } 236403adfb6SMatthew G Knepley 237a6dfd86eSKarl Rupp bool STREAM_Copy_verify_double(double *a, double *b, size_t len) 238a6dfd86eSKarl Rupp { 239caccb7e3SMatthew G Knepley size_t idx; 240caccb7e3SMatthew G Knepley bool bDifferent = false; 241caccb7e3SMatthew G Knepley 242caccb7e3SMatthew G Knepley for (idx = 0; idx < len && !bDifferent; idx++) { 243caccb7e3SMatthew G Knepley double expectedResult = a[idx]; 244caccb7e3SMatthew G Knepley double diffResultExpected = (b[idx] - expectedResult); 24519816777SMark double relErrorULPS = (fabsf(diffResultExpected)/fabsf(expectedResult))/dbl_eps; 246caccb7e3SMatthew G Knepley /* element-wise relative error determination */ 247caccb7e3SMatthew G Knepley bDifferent = (relErrorULPS > 2.); 248caccb7e3SMatthew G Knepley } 249caccb7e3SMatthew G Knepley 250caccb7e3SMatthew G Knepley return bDifferent; 251caccb7e3SMatthew G Knepley } 252caccb7e3SMatthew G Knepley 253a6dfd86eSKarl Rupp bool STREAM_Scale_verify(float *a, float *b, float scale, size_t len) 254a6dfd86eSKarl Rupp { 255403adfb6SMatthew G Knepley size_t idx; 256403adfb6SMatthew G Knepley bool bDifferent = false; 257403adfb6SMatthew G Knepley 258403adfb6SMatthew G Knepley for (idx = 0; idx < len && !bDifferent; idx++) { 259403adfb6SMatthew G Knepley float expectedResult = scale*a[idx]; 260403adfb6SMatthew G Knepley float diffResultExpected = (b[idx] - expectedResult); 261403adfb6SMatthew G Knepley float relErrorULPS = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps; 262403adfb6SMatthew G Knepley /* element-wise relative error determination */ 263403adfb6SMatthew G Knepley bDifferent = (relErrorULPS > 2.f); 264403adfb6SMatthew G Knepley } 265403adfb6SMatthew G Knepley 266403adfb6SMatthew G Knepley return bDifferent; 267403adfb6SMatthew G Knepley } 268403adfb6SMatthew G Knepley 269a6dfd86eSKarl Rupp bool STREAM_Scale_verify_double(double *a, double *b, double scale, size_t len) 270a6dfd86eSKarl Rupp { 271caccb7e3SMatthew G Knepley size_t idx; 272caccb7e3SMatthew G Knepley bool bDifferent = false; 273caccb7e3SMatthew G Knepley 274caccb7e3SMatthew G Knepley for (idx = 0; idx < len && !bDifferent; idx++) { 275caccb7e3SMatthew G Knepley double expectedResult = scale*a[idx]; 276caccb7e3SMatthew G Knepley double diffResultExpected = (b[idx] - expectedResult); 277caccb7e3SMatthew G Knepley double relErrorULPS = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps; 278caccb7e3SMatthew G Knepley /* element-wise relative error determination */ 279caccb7e3SMatthew G Knepley bDifferent = (relErrorULPS > 2.); 280caccb7e3SMatthew G Knepley } 281caccb7e3SMatthew G Knepley 282caccb7e3SMatthew G Knepley return bDifferent; 283caccb7e3SMatthew G Knepley } 284caccb7e3SMatthew G Knepley 285a6dfd86eSKarl Rupp bool STREAM_Add_verify(float *a, float *b, float *c, size_t len) 286a6dfd86eSKarl Rupp { 287403adfb6SMatthew G Knepley size_t idx; 288403adfb6SMatthew G Knepley bool bDifferent = false; 289403adfb6SMatthew G Knepley 290403adfb6SMatthew G Knepley for (idx = 0; idx < len && !bDifferent; idx++) { 291403adfb6SMatthew G Knepley float expectedResult = a[idx] + b[idx]; 292403adfb6SMatthew G Knepley float diffResultExpected = (c[idx] - expectedResult); 293403adfb6SMatthew G Knepley float relErrorULPS = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps; 294403adfb6SMatthew G Knepley /* element-wise relative error determination */ 295403adfb6SMatthew G Knepley bDifferent = (relErrorULPS > 2.f); 296403adfb6SMatthew G Knepley } 297403adfb6SMatthew G Knepley 298403adfb6SMatthew G Knepley return bDifferent; 299403adfb6SMatthew G Knepley } 300403adfb6SMatthew G Knepley 301a6dfd86eSKarl Rupp bool STREAM_Add_verify_double(double *a, double *b, double *c, size_t len) 302a6dfd86eSKarl Rupp { 303caccb7e3SMatthew G Knepley size_t idx; 304caccb7e3SMatthew G Knepley bool bDifferent = false; 305caccb7e3SMatthew G Knepley 306caccb7e3SMatthew G Knepley for (idx = 0; idx < len && !bDifferent; idx++) { 307caccb7e3SMatthew G Knepley double expectedResult = a[idx] + b[idx]; 308caccb7e3SMatthew G Knepley double diffResultExpected = (c[idx] - expectedResult); 309caccb7e3SMatthew G Knepley double relErrorULPS = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps; 310caccb7e3SMatthew G Knepley /* element-wise relative error determination */ 311caccb7e3SMatthew G Knepley bDifferent = (relErrorULPS > 2.); 312caccb7e3SMatthew G Knepley } 313caccb7e3SMatthew G Knepley 314caccb7e3SMatthew G Knepley return bDifferent; 315caccb7e3SMatthew G Knepley } 316caccb7e3SMatthew G Knepley 317a6dfd86eSKarl Rupp bool STREAM_Triad_verify(float *a, float *b, float *c, float scalar, size_t len) 318a6dfd86eSKarl Rupp { 319403adfb6SMatthew G Knepley size_t idx; 320403adfb6SMatthew G Knepley bool bDifferent = false; 321403adfb6SMatthew G Knepley 322403adfb6SMatthew G Knepley for (idx = 0; idx < len && !bDifferent; idx++) { 323403adfb6SMatthew G Knepley float expectedResult = a[idx] + scalar*b[idx]; 324403adfb6SMatthew G Knepley float diffResultExpected = (c[idx] - expectedResult); 325403adfb6SMatthew G Knepley float relErrorULPS = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps; 326403adfb6SMatthew G Knepley /* element-wise relative error determination */ 327403adfb6SMatthew G Knepley bDifferent = (relErrorULPS > 3.f); 328403adfb6SMatthew G Knepley } 329403adfb6SMatthew G Knepley 330403adfb6SMatthew G Knepley return bDifferent; 331403adfb6SMatthew G Knepley } 332403adfb6SMatthew G Knepley 333a6dfd86eSKarl Rupp bool STREAM_Triad_verify_double(double *a, double *b, double *c, double scalar, size_t len) 334a6dfd86eSKarl Rupp { 335caccb7e3SMatthew G Knepley size_t idx; 336caccb7e3SMatthew G Knepley bool bDifferent = false; 337caccb7e3SMatthew G Knepley 338caccb7e3SMatthew G Knepley for (idx = 0; idx < len && !bDifferent; idx++) { 339caccb7e3SMatthew G Knepley double expectedResult = a[idx] + scalar*b[idx]; 340caccb7e3SMatthew G Knepley double diffResultExpected = (c[idx] - expectedResult); 341caccb7e3SMatthew G Knepley double relErrorULPS = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps; 342caccb7e3SMatthew G Knepley /* element-wise relative error determination */ 343caccb7e3SMatthew G Knepley bDifferent = (relErrorULPS > 3.); 344caccb7e3SMatthew G Knepley } 345caccb7e3SMatthew G Knepley 346caccb7e3SMatthew G Knepley return bDifferent; 347caccb7e3SMatthew G Knepley } 348caccb7e3SMatthew G Knepley 349403adfb6SMatthew G Knepley /* forward declarations */ 350caccb7e3SMatthew G Knepley PetscErrorCode setupStream(PetscInt device, PetscBool runDouble, PetscBool cpuTiming); 351403adfb6SMatthew G Knepley PetscErrorCode runStream(const PetscInt iNumThreadsPerBlock, PetscBool bDontUseGPUTiming); 352caccb7e3SMatthew G Knepley PetscErrorCode runStreamDouble(const PetscInt iNumThreadsPerBlock, PetscBool bDontUseGPUTiming); 35319816777SMark PetscErrorCode printResultsReadable(float times[][NTIMES], size_t); 354403adfb6SMatthew G Knepley 355403adfb6SMatthew G Knepley int main(int argc, char *argv[]) 356403adfb6SMatthew G Knepley { 357403adfb6SMatthew G Knepley PetscInt device = 0; 35819816777SMark PetscBool runDouble = PETSC_TRUE; 35919816777SMark const PetscBool cpuTiming = PETSC_TRUE; // must be true 360403adfb6SMatthew G Knepley PetscErrorCode ierr; 361403adfb6SMatthew G Knepley 3629566063dSJacob Faibussowitsch PetscCallCUDA(cudaSetDeviceFlags(cudaDeviceBlockingSync)); 36319816777SMark 3649566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 365403adfb6SMatthew G Knepley 366d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, "", "STREAM Benchmark Options", "STREAM"); 3679566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-device", "Specify the CUDA device to be used", "STREAM", device, &device, NULL,0)); 3689566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-double", "Also run double precision tests", "STREAM", runDouble, &runDouble, NULL)); 369d0609cedSBarry Smith PetscOptionsEnd(); 370403adfb6SMatthew G Knepley 371caccb7e3SMatthew G Knepley ierr = setupStream(device, runDouble, cpuTiming); 37219816777SMark if (ierr) { 3739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[streamBenchmark] - results:\t%s\n\n", (ierr == 0) ? "PASSES" : "FAILED")); 374403adfb6SMatthew G Knepley } 3759566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 376b122ec5aSJacob Faibussowitsch return 0; 377403adfb6SMatthew G Knepley } 378403adfb6SMatthew G Knepley 379403adfb6SMatthew G Knepley /////////////////////////////////////////////////////////////////////////////// 380403adfb6SMatthew G Knepley //Run the appropriate tests 381403adfb6SMatthew G Knepley /////////////////////////////////////////////////////////////////////////////// 382caccb7e3SMatthew G Knepley PetscErrorCode setupStream(PetscInt deviceNum, PetscBool runDouble, PetscBool cpuTiming) 383403adfb6SMatthew G Knepley { 384403adfb6SMatthew G Knepley PetscInt iNumThreadsPerBlock = 128; 385403adfb6SMatthew G Knepley 386403adfb6SMatthew G Knepley PetscFunctionBegin; 387403adfb6SMatthew G Knepley // Check device 388403adfb6SMatthew G Knepley { 389403adfb6SMatthew G Knepley int deviceCount; 390403adfb6SMatthew G Knepley 391403adfb6SMatthew G Knepley cudaGetDeviceCount(&deviceCount); 392403adfb6SMatthew G Knepley if (deviceCount == 0) { 3939566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "!!!!!No devices found!!!!!\n")); 394403adfb6SMatthew G Knepley return -1000; 395403adfb6SMatthew G Knepley } 396403adfb6SMatthew G Knepley 397403adfb6SMatthew G Knepley if (deviceNum >= deviceCount || deviceNum < 0) { 3989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n!!!!!Invalid GPU number %d given hence default gpu %d will be used !!!!!\n", deviceNum, 0)); 399403adfb6SMatthew G Knepley deviceNum = 0; 400403adfb6SMatthew G Knepley } 401403adfb6SMatthew G Knepley } 402403adfb6SMatthew G Knepley 403403adfb6SMatthew G Knepley cudaSetDevice(deviceNum); 4049566063dSJacob Faibussowitsch // PetscCall(PetscPrintf(PETSC_COMM_SELF, "Running on...\n\n")); 405403adfb6SMatthew G Knepley cudaDeviceProp deviceProp; 40619816777SMark if (cudaGetDeviceProperties(&deviceProp, deviceNum) != cudaSuccess) { 4079566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Unable to determine device %d properties, exiting\n")); 408403adfb6SMatthew G Knepley return -1; 409403adfb6SMatthew G Knepley } 410403adfb6SMatthew G Knepley 411caccb7e3SMatthew G Knepley if (runDouble && deviceProp.major == 1 && deviceProp.minor < 3) { 4129566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Unable to run double-precision STREAM benchmark on a compute capability GPU less than 1.3\n")); 413caccb7e3SMatthew G Knepley return -1; 414caccb7e3SMatthew G Knepley } 4156f2b61bcSKarl Rupp if (deviceProp.major == 2 && deviceProp.minor == 1) iNumThreadsPerBlock = 192; /* GF104 architecture / 48 CUDA Cores per MP */ 4166f2b61bcSKarl Rupp else iNumThreadsPerBlock = 128; /* GF100 architecture / 32 CUDA Cores per MP */ 417403adfb6SMatthew G Knepley 418*1baa6e33SBarry Smith if (runDouble) PetscCall(runStreamDouble(iNumThreadsPerBlock, cpuTiming)); 419*1baa6e33SBarry Smith else PetscCall(runStream(iNumThreadsPerBlock, cpuTiming)); 420403adfb6SMatthew G Knepley PetscFunctionReturn(0); 421403adfb6SMatthew G Knepley } 422403adfb6SMatthew G Knepley 423403adfb6SMatthew G Knepley /////////////////////////////////////////////////////////////////////////// 424403adfb6SMatthew G Knepley // runStream 425403adfb6SMatthew G Knepley /////////////////////////////////////////////////////////////////////////// 426403adfb6SMatthew G Knepley PetscErrorCode runStream(const PetscInt iNumThreadsPerBlock, PetscBool bDontUseGPUTiming) 427403adfb6SMatthew G Knepley { 428403adfb6SMatthew G Knepley float *d_a, *d_b, *d_c; 429403adfb6SMatthew G Knepley int k; 430caccb7e3SMatthew G Knepley float times[8][NTIMES]; 431403adfb6SMatthew G Knepley float scalar; 432403adfb6SMatthew G Knepley 433403adfb6SMatthew G Knepley PetscFunctionBegin; 434403adfb6SMatthew G Knepley /* Allocate memory on device */ 4359566063dSJacob Faibussowitsch PetscCallCUDA(cudaMalloc((void**)&d_a, sizeof(float)*N)); 4369566063dSJacob Faibussowitsch PetscCallCUDA(cudaMalloc((void**)&d_b, sizeof(float)*N)); 4379566063dSJacob Faibussowitsch PetscCallCUDA(cudaMalloc((void**)&d_c, sizeof(float)*N)); 438403adfb6SMatthew G Knepley 439403adfb6SMatthew G Knepley /* Compute execution configuration */ 440403adfb6SMatthew G Knepley 441403adfb6SMatthew G Knepley dim3 dimBlock(iNumThreadsPerBlock); /* (iNumThreadsPerBlock,1,1) */ 442403adfb6SMatthew G Knepley dim3 dimGrid(N/dimBlock.x); /* (N/dimBlock.x,1,1) */ 443403adfb6SMatthew G Knepley if (N % dimBlock.x != 0) dimGrid.x+=1; 444403adfb6SMatthew G Knepley 445403adfb6SMatthew G Knepley /* Initialize memory on the device */ 446403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_a, 2.f, N); 447403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N); 448403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N); 449403adfb6SMatthew G Knepley 450403adfb6SMatthew G Knepley /* --- MAIN LOOP --- repeat test cases NTIMES times --- */ 451403adfb6SMatthew G Knepley PetscLogDouble cpuTimer = 0.0; 452403adfb6SMatthew G Knepley 453403adfb6SMatthew G Knepley scalar=3.0f; 454403adfb6SMatthew G Knepley for (k = 0; k < NTIMES; ++k) { 4558563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 456403adfb6SMatthew G Knepley STREAM_Copy<<<dimGrid,dimBlock>>>(d_a, d_c, N); 45719816777SMark cudaStreamSynchronize(NULL); 4589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 4598563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 46019816777SMark if (bDontUseGPUTiming) times[0][k] = cpuTimer*1.e3; // millisec 461403adfb6SMatthew G Knepley 462403adfb6SMatthew G Knepley cpuTimer = 0.0; 4638563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 464403adfb6SMatthew G Knepley STREAM_Copy_Optimized<<<dimGrid,dimBlock>>>(d_a, d_c, N); 46519816777SMark cudaStreamSynchronize(NULL); 4669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 467df3898eeSBarry Smith //get the total elapsed time in ms 4688563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 46919816777SMark if (bDontUseGPUTiming) times[1][k] = cpuTimer*1.e3; 470403adfb6SMatthew G Knepley 471403adfb6SMatthew G Knepley cpuTimer = 0.0; 4728563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 473403adfb6SMatthew G Knepley STREAM_Scale<<<dimGrid,dimBlock>>>(d_b, d_c, scalar, N); 47419816777SMark cudaStreamSynchronize(NULL); 4759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 476df3898eeSBarry Smith //get the total elapsed time in ms 4778563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 47819816777SMark if (bDontUseGPUTiming) times[2][k] = cpuTimer*1.e3; 479403adfb6SMatthew G Knepley 480403adfb6SMatthew G Knepley cpuTimer = 0.0; 4818563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 482caccb7e3SMatthew G Knepley STREAM_Scale_Optimized<<<dimGrid,dimBlock>>>(d_b, d_c, scalar, N); 48319816777SMark cudaStreamSynchronize(NULL); 4849566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 485df3898eeSBarry Smith //get the total elapsed time in ms 4868563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 48719816777SMark if (bDontUseGPUTiming) times[3][k] = cpuTimer*1.e3; 488403adfb6SMatthew G Knepley 489403adfb6SMatthew G Knepley cpuTimer = 0.0; 4908563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 4919566063dSJacob Faibussowitsch // PetscCallCUDA(cudaEventRecord(start, 0)); 492caccb7e3SMatthew G Knepley STREAM_Add<<<dimGrid,dimBlock>>>(d_a, d_b, d_c, N); 49319816777SMark cudaStreamSynchronize(NULL); 4949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 4959566063dSJacob Faibussowitsch PetscCallCUDA(cudaEventRecord(stop, 0)); 4969566063dSJacob Faibussowitsch // PetscCallCUDA(cudaEventSynchronize(stop)); 497df3898eeSBarry Smith //get the total elapsed time in ms 4988563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 49919816777SMark if (bDontUseGPUTiming) times[4][k] = cpuTimer*1.e3; 5006f2b61bcSKarl Rupp else { 5019566063dSJacob Faibussowitsch // PetscCallCUDA(cudaEventElapsedTime(×[4][k], start, stop)); 502403adfb6SMatthew G Knepley } 503403adfb6SMatthew G Knepley 504caccb7e3SMatthew G Knepley cpuTimer = 0.0; 5058563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 506caccb7e3SMatthew G Knepley STREAM_Add_Optimized<<<dimGrid,dimBlock>>>(d_a, d_b, d_c, N); 50719816777SMark cudaStreamSynchronize(NULL); 5089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 509df3898eeSBarry Smith //get the total elapsed time in ms 5108563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 51119816777SMark if (bDontUseGPUTiming) times[5][k] = cpuTimer*1.e3; 512caccb7e3SMatthew G Knepley 513caccb7e3SMatthew G Knepley cpuTimer = 0.0; 5148563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 515caccb7e3SMatthew G Knepley STREAM_Triad<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar, N); 51619816777SMark cudaStreamSynchronize(NULL); 5179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 518df3898eeSBarry Smith //get the total elapsed time in ms 5198563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 52019816777SMark if (bDontUseGPUTiming) times[6][k] = cpuTimer*1.e3; 521caccb7e3SMatthew G Knepley 522caccb7e3SMatthew G Knepley cpuTimer = 0.0; 5238563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 524caccb7e3SMatthew G Knepley STREAM_Triad_Optimized<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar, N); 52519816777SMark cudaStreamSynchronize(NULL); 5269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 527df3898eeSBarry Smith //get the total elapsed time in ms 5288563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 52919816777SMark if (bDontUseGPUTiming) times[7][k] = cpuTimer*1.e3; 530caccb7e3SMatthew G Knepley } 531caccb7e3SMatthew G Knepley 53219816777SMark if (1) { /* verify kernels */ 533403adfb6SMatthew G Knepley float *h_a, *h_b, *h_c; 534403adfb6SMatthew G Knepley bool errorSTREAMkernel = true; 535403adfb6SMatthew G Knepley 536403adfb6SMatthew G Knepley if ((h_a = (float*)calloc(N, sizeof(float))) == (float*)NULL) { 537403adfb6SMatthew G Knepley printf("Unable to allocate array h_a, exiting ...\n"); 538403adfb6SMatthew G Knepley exit(1); 539403adfb6SMatthew G Knepley } 540403adfb6SMatthew G Knepley if ((h_b = (float*)calloc(N, sizeof(float))) == (float*)NULL) { 541403adfb6SMatthew G Knepley printf("Unable to allocate array h_b, exiting ...\n"); 542403adfb6SMatthew G Knepley exit(1); 543403adfb6SMatthew G Knepley } 544403adfb6SMatthew G Knepley 545403adfb6SMatthew G Knepley if ((h_c = (float*)calloc(N, sizeof(float))) == (float*)NULL) { 546403adfb6SMatthew G Knepley printf("Unalbe to allocate array h_c, exiting ...\n"); 547403adfb6SMatthew G Knepley exit(1); 548403adfb6SMatthew G Knepley } 549403adfb6SMatthew G Knepley 550403adfb6SMatthew G Knepley /* 551403adfb6SMatthew G Knepley * perform kernel, copy device memory into host memory and verify each 552403adfb6SMatthew G Knepley * device kernel output 553403adfb6SMatthew G Knepley */ 554403adfb6SMatthew G Knepley 555403adfb6SMatthew G Knepley /* Initialize memory on the device */ 556403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_a, 2.f, N); 557403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N); 558403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N); 559403adfb6SMatthew G Knepley 560403adfb6SMatthew G Knepley STREAM_Copy<<<dimGrid,dimBlock>>>(d_a, d_c, N); 5619566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_a, d_a, sizeof(float) * N, cudaMemcpyDeviceToHost)); 5629566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_c, d_c, sizeof(float) * N, cudaMemcpyDeviceToHost)); 563403adfb6SMatthew G Knepley errorSTREAMkernel = STREAM_Copy_verify(h_a, h_c, N); 564403adfb6SMatthew G Knepley if (errorSTREAMkernel) { 5659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy:\t\tError detected in device STREAM_Copy, exiting\n")); 566403adfb6SMatthew G Knepley exit(-2000); 567403adfb6SMatthew G Knepley } 568403adfb6SMatthew G Knepley 569403adfb6SMatthew G Knepley /* Initialize memory on the device */ 570403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_a, 2.f, N); 571403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N); 572403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N); 573403adfb6SMatthew G Knepley 574403adfb6SMatthew G Knepley STREAM_Copy_Optimized<<<dimGrid,dimBlock>>>(d_a, d_c, N); 5759566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_a, d_a, sizeof(float) * N, cudaMemcpyDeviceToHost)); 5769566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_c, d_c, sizeof(float) * N, cudaMemcpyDeviceToHost)); 577403adfb6SMatthew G Knepley errorSTREAMkernel = STREAM_Copy_verify(h_a, h_c, N); 578403adfb6SMatthew G Knepley if (errorSTREAMkernel) { 5799566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy_Optimized:\tError detected in device STREAM_Copy_Optimized, exiting\n")); 580403adfb6SMatthew G Knepley exit(-3000); 581403adfb6SMatthew G Knepley } 582403adfb6SMatthew G Knepley 583403adfb6SMatthew G Knepley /* Initialize memory on the device */ 58419816777SMark set_array<<<dimGrid,dimBlock>>>(d_a, 2.f, N); 585403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N); 586403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N); 587403adfb6SMatthew G Knepley 588403adfb6SMatthew G Knepley STREAM_Scale<<<dimGrid,dimBlock>>>(d_b, d_c, scalar, N); 5899566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_b, d_b, sizeof(float) * N, cudaMemcpyDeviceToHost)); 5909566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_c, d_c, sizeof(float) * N, cudaMemcpyDeviceToHost)); 591403adfb6SMatthew G Knepley errorSTREAMkernel = STREAM_Scale_verify(h_b, h_c, scalar, N); 592403adfb6SMatthew G Knepley if (errorSTREAMkernel) { 5939566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " device STREAM_Scale:\t\tError detected in device STREAM_Scale, exiting\n")); 594403adfb6SMatthew G Knepley exit(-4000); 595403adfb6SMatthew G Knepley } 596403adfb6SMatthew G Knepley 597403adfb6SMatthew G Knepley /* Initialize memory on the device */ 598403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_a, 2.f, N); 599403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N); 600403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N); 601403adfb6SMatthew G Knepley 602403adfb6SMatthew G Knepley STREAM_Add<<<dimGrid,dimBlock>>>(d_a, d_b, d_c, N); 6039566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_a, d_a, sizeof(float) * N, cudaMemcpyDeviceToHost)); 6049566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_b, d_b, sizeof(float) * N, cudaMemcpyDeviceToHost)); 6059566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_c, d_c, sizeof(float) * N, cudaMemcpyDeviceToHost)); 606403adfb6SMatthew G Knepley errorSTREAMkernel = STREAM_Add_verify(h_a, h_b, h_c, N); 607403adfb6SMatthew G Knepley if (errorSTREAMkernel) { 6089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " device STREAM_Add:\t\tError detected in device STREAM_Add, exiting\n")); 609403adfb6SMatthew G Knepley exit(-5000); 610403adfb6SMatthew G Knepley } 611403adfb6SMatthew G Knepley 612403adfb6SMatthew G Knepley /* Initialize memory on the device */ 613403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_a, 2.f, N); 614403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N); 615403adfb6SMatthew G Knepley set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N); 616403adfb6SMatthew G Knepley 617403adfb6SMatthew G Knepley STREAM_Triad<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar, N); 6189566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_a, d_a, sizeof(float) * N, cudaMemcpyDeviceToHost)); 6199566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_b, d_b, sizeof(float) * N, cudaMemcpyDeviceToHost)); 6209566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_c, d_c, sizeof(float) * N, cudaMemcpyDeviceToHost)); 621403adfb6SMatthew G Knepley errorSTREAMkernel = STREAM_Triad_verify(h_b, h_c, h_a, scalar, N); 622403adfb6SMatthew G Knepley if (errorSTREAMkernel) { 6239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " device STREAM_Triad:\t\tError detected in device STREAM_Triad, exiting\n")); 624403adfb6SMatthew G Knepley exit(-6000); 625403adfb6SMatthew G Knepley } 626403adfb6SMatthew G Knepley 62719816777SMark free(h_a); 62819816777SMark free(h_b); 62919816777SMark free(h_c); 63019816777SMark } 631403adfb6SMatthew G Knepley /* continue from here */ 63219816777SMark printResultsReadable(times, sizeof(float)); 633403adfb6SMatthew G Knepley 634403adfb6SMatthew G Knepley /* Free memory on device */ 6359566063dSJacob Faibussowitsch PetscCallCUDA(cudaFree(d_a)); 6369566063dSJacob Faibussowitsch PetscCallCUDA(cudaFree(d_b)); 6379566063dSJacob Faibussowitsch PetscCallCUDA(cudaFree(d_c)); 63819816777SMark 639403adfb6SMatthew G Knepley PetscFunctionReturn(0); 640403adfb6SMatthew G Knepley } 641403adfb6SMatthew G Knepley 642caccb7e3SMatthew G Knepley PetscErrorCode runStreamDouble(const PetscInt iNumThreadsPerBlock, PetscBool bDontUseGPUTiming) 643caccb7e3SMatthew G Knepley { 644caccb7e3SMatthew G Knepley double *d_a, *d_b, *d_c; 645caccb7e3SMatthew G Knepley int k; 646caccb7e3SMatthew G Knepley float times[8][NTIMES]; 647caccb7e3SMatthew G Knepley double scalar; 648caccb7e3SMatthew G Knepley 649caccb7e3SMatthew G Knepley PetscFunctionBegin; 650caccb7e3SMatthew G Knepley /* Allocate memory on device */ 6519566063dSJacob Faibussowitsch PetscCallCUDA(cudaMalloc((void**)&d_a, sizeof(double)*N)); 6529566063dSJacob Faibussowitsch PetscCallCUDA(cudaMalloc((void**)&d_b, sizeof(double)*N)); 6539566063dSJacob Faibussowitsch PetscCallCUDA(cudaMalloc((void**)&d_c, sizeof(double)*N)); 654caccb7e3SMatthew G Knepley 655caccb7e3SMatthew G Knepley /* Compute execution configuration */ 656caccb7e3SMatthew G Knepley 657caccb7e3SMatthew G Knepley dim3 dimBlock(iNumThreadsPerBlock); /* (iNumThreadsPerBlock,1,1) */ 658caccb7e3SMatthew G Knepley dim3 dimGrid(N/dimBlock.x); /* (N/dimBlock.x,1,1) */ 659caccb7e3SMatthew G Knepley if (N % dimBlock.x != 0) dimGrid.x+=1; 660caccb7e3SMatthew G Knepley 661caccb7e3SMatthew G Knepley /* Initialize memory on the device */ 662caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_a, 2., N); 663caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N); 664caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N); 665caccb7e3SMatthew G Knepley 666caccb7e3SMatthew G Knepley /* --- MAIN LOOP --- repeat test cases NTIMES times --- */ 667caccb7e3SMatthew G Knepley PetscLogDouble cpuTimer = 0.0; 668caccb7e3SMatthew G Knepley 669caccb7e3SMatthew G Knepley scalar=3.0; 670caccb7e3SMatthew G Knepley for (k = 0; k < NTIMES; ++k) { 6718563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 672caccb7e3SMatthew G Knepley STREAM_Copy_double<<<dimGrid,dimBlock>>>(d_a, d_c, N); 67319816777SMark cudaStreamSynchronize(NULL); 6749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 675df3898eeSBarry Smith //get the total elapsed time in ms 676caccb7e3SMatthew G Knepley if (bDontUseGPUTiming) { 6778563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 67819816777SMark times[0][k] = cpuTimer*1.e3; 679caccb7e3SMatthew G Knepley } 680caccb7e3SMatthew G Knepley 681caccb7e3SMatthew G Knepley cpuTimer = 0.0; 6828563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 683caccb7e3SMatthew G Knepley STREAM_Copy_Optimized_double<<<dimGrid,dimBlock>>>(d_a, d_c, N); 68419816777SMark cudaStreamSynchronize(NULL); 6859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 686df3898eeSBarry Smith //get the total elapsed time in ms 687caccb7e3SMatthew G Knepley if (bDontUseGPUTiming) { 6888563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 68919816777SMark times[1][k] = cpuTimer*1.e3; 690caccb7e3SMatthew G Knepley } 691caccb7e3SMatthew G Knepley 692caccb7e3SMatthew G Knepley cpuTimer = 0.0; 6938563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 694caccb7e3SMatthew G Knepley STREAM_Scale_double<<<dimGrid,dimBlock>>>(d_b, d_c, scalar, N); 69519816777SMark cudaStreamSynchronize(NULL); 6969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 697df3898eeSBarry Smith //get the total elapsed time in ms 6988563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 69919816777SMark if (bDontUseGPUTiming) times[2][k] = cpuTimer*1.e3; 700caccb7e3SMatthew G Knepley 701caccb7e3SMatthew G Knepley cpuTimer = 0.0; 7028563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 703caccb7e3SMatthew G Knepley STREAM_Scale_Optimized_double<<<dimGrid,dimBlock>>>(d_b, d_c, scalar, N); 70419816777SMark cudaStreamSynchronize(NULL); 7059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 706df3898eeSBarry Smith //get the total elapsed time in ms 7078563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 70819816777SMark if (bDontUseGPUTiming) times[3][k] = cpuTimer*1.e3; 709caccb7e3SMatthew G Knepley 710caccb7e3SMatthew G Knepley cpuTimer = 0.0; 7118563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 712caccb7e3SMatthew G Knepley STREAM_Add_double<<<dimGrid,dimBlock>>>(d_a, d_b, d_c, N); 71319816777SMark cudaStreamSynchronize(NULL); 7149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 715df3898eeSBarry Smith //get the total elapsed time in ms 7168563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 71719816777SMark if (bDontUseGPUTiming) times[4][k] = cpuTimer*1.e3; 718caccb7e3SMatthew G Knepley 719caccb7e3SMatthew G Knepley cpuTimer = 0.0; 7208563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 721caccb7e3SMatthew G Knepley STREAM_Add_Optimized_double<<<dimGrid,dimBlock>>>(d_a, d_b, d_c, N); 72219816777SMark cudaStreamSynchronize(NULL); 7239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 724df3898eeSBarry Smith //get the total elapsed time in ms 7258563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 72619816777SMark if (bDontUseGPUTiming) times[5][k] = cpuTimer*1.e3; 727caccb7e3SMatthew G Knepley 728caccb7e3SMatthew G Knepley cpuTimer = 0.0; 7298563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 730caccb7e3SMatthew G Knepley STREAM_Triad_double<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar, N); 73119816777SMark cudaStreamSynchronize(NULL); 7329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 733df3898eeSBarry Smith //get the total elapsed time in ms 7348563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 73519816777SMark if (bDontUseGPUTiming) times[6][k] = cpuTimer*1.e3; 736caccb7e3SMatthew G Knepley 737caccb7e3SMatthew G Knepley cpuTimer = 0.0; 7388563dfccSBarry Smith PetscTimeSubtract(&cpuTimer); 739caccb7e3SMatthew G Knepley STREAM_Triad_Optimized_double<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar, N); 74019816777SMark cudaStreamSynchronize(NULL); 7419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(MPI_COMM_WORLD)); 742df3898eeSBarry Smith //get the total elapsed time in ms 7438563dfccSBarry Smith PetscTimeAdd(&cpuTimer); 74419816777SMark if (bDontUseGPUTiming) times[7][k] = cpuTimer*1.e3; 745caccb7e3SMatthew G Knepley } 746caccb7e3SMatthew G Knepley 74719816777SMark if (1) { /* verify kernels */ 748caccb7e3SMatthew G Knepley double *h_a, *h_b, *h_c; 749caccb7e3SMatthew G Knepley bool errorSTREAMkernel = true; 750caccb7e3SMatthew G Knepley 751caccb7e3SMatthew G Knepley if ((h_a = (double*)calloc(N, sizeof(double))) == (double*)NULL) { 752caccb7e3SMatthew G Knepley printf("Unable to allocate array h_a, exiting ...\n"); 753caccb7e3SMatthew G Knepley exit(1); 754caccb7e3SMatthew G Knepley } 755caccb7e3SMatthew G Knepley if ((h_b = (double*)calloc(N, sizeof(double))) == (double*)NULL) { 756caccb7e3SMatthew G Knepley printf("Unable to allocate array h_b, exiting ...\n"); 757caccb7e3SMatthew G Knepley exit(1); 758caccb7e3SMatthew G Knepley } 759caccb7e3SMatthew G Knepley 760caccb7e3SMatthew G Knepley if ((h_c = (double*)calloc(N, sizeof(double))) == (double*)NULL) { 761caccb7e3SMatthew G Knepley printf("Unalbe to allocate array h_c, exiting ...\n"); 762caccb7e3SMatthew G Knepley exit(1); 763caccb7e3SMatthew G Knepley } 764caccb7e3SMatthew G Knepley 765caccb7e3SMatthew G Knepley /* 766caccb7e3SMatthew G Knepley * perform kernel, copy device memory into host memory and verify each 767caccb7e3SMatthew G Knepley * device kernel output 768caccb7e3SMatthew G Knepley */ 769caccb7e3SMatthew G Knepley 770caccb7e3SMatthew G Knepley /* Initialize memory on the device */ 771caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_a, 2., N); 772caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N); 773caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N); 774caccb7e3SMatthew G Knepley 775caccb7e3SMatthew G Knepley STREAM_Copy_double<<<dimGrid,dimBlock>>>(d_a, d_c, N); 7769566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_a, d_a, sizeof(double) * N, cudaMemcpyDeviceToHost)); 7779566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_c, d_c, sizeof(double) * N, cudaMemcpyDeviceToHost)); 778caccb7e3SMatthew G Knepley errorSTREAMkernel = STREAM_Copy_verify_double(h_a, h_c, N); 779caccb7e3SMatthew G Knepley if (errorSTREAMkernel) { 7809566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy:\t\tError detected in device STREAM_Copy, exiting\n")); 781caccb7e3SMatthew G Knepley exit(-2000); 782caccb7e3SMatthew G Knepley } 783caccb7e3SMatthew G Knepley 784caccb7e3SMatthew G Knepley /* Initialize memory on the device */ 785caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_a, 2., N); 786caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N); 787caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N); 788caccb7e3SMatthew G Knepley 789caccb7e3SMatthew G Knepley STREAM_Copy_Optimized_double<<<dimGrid,dimBlock>>>(d_a, d_c, N); 7909566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_a, d_a, sizeof(double) * N, cudaMemcpyDeviceToHost)); 7919566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_c, d_c, sizeof(double) * N, cudaMemcpyDeviceToHost)); 792caccb7e3SMatthew G Knepley errorSTREAMkernel = STREAM_Copy_verify_double(h_a, h_c, N); 793caccb7e3SMatthew G Knepley if (errorSTREAMkernel) { 7949566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy_Optimized:\tError detected in device STREAM_Copy_Optimized, exiting\n")); 795caccb7e3SMatthew G Knepley exit(-3000); 796caccb7e3SMatthew G Knepley } 797caccb7e3SMatthew G Knepley 798caccb7e3SMatthew G Knepley /* Initialize memory on the device */ 799caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N); 800caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N); 801caccb7e3SMatthew G Knepley 802caccb7e3SMatthew G Knepley STREAM_Scale_double<<<dimGrid,dimBlock>>>(d_b, d_c, scalar, N); 8039566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_b, d_b, sizeof(double) * N, cudaMemcpyDeviceToHost)); 8049566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_c, d_c, sizeof(double) * N, cudaMemcpyDeviceToHost)); 805caccb7e3SMatthew G Knepley errorSTREAMkernel = STREAM_Scale_verify_double(h_b, h_c, scalar, N); 806caccb7e3SMatthew G Knepley if (errorSTREAMkernel) { 8079566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " device STREAM_Scale:\t\tError detected in device STREAM_Scale, exiting\n")); 808caccb7e3SMatthew G Knepley exit(-4000); 809caccb7e3SMatthew G Knepley } 810caccb7e3SMatthew G Knepley 811caccb7e3SMatthew G Knepley /* Initialize memory on the device */ 812caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_a, 2., N); 813caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N); 814caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N); 815caccb7e3SMatthew G Knepley 816caccb7e3SMatthew G Knepley STREAM_Add_double<<<dimGrid,dimBlock>>>(d_a, d_b, d_c, N); 8179566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_a, d_a, sizeof(double) * N, cudaMemcpyDeviceToHost)); 8189566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_b, d_b, sizeof(double) * N, cudaMemcpyDeviceToHost)); 8199566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_c, d_c, sizeof(double) * N, cudaMemcpyDeviceToHost)); 820caccb7e3SMatthew G Knepley errorSTREAMkernel = STREAM_Add_verify_double(h_a, h_b, h_c, N); 821caccb7e3SMatthew G Knepley if (errorSTREAMkernel) { 8229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " device STREAM_Add:\t\tError detected in device STREAM_Add, exiting\n")); 823caccb7e3SMatthew G Knepley exit(-5000); 824caccb7e3SMatthew G Knepley } 825caccb7e3SMatthew G Knepley 826caccb7e3SMatthew G Knepley /* Initialize memory on the device */ 827caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_a, 2., N); 828caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N); 829caccb7e3SMatthew G Knepley set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N); 830caccb7e3SMatthew G Knepley 831caccb7e3SMatthew G Knepley STREAM_Triad_double<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar, N); 8329566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_a, d_a, sizeof(double) * N, cudaMemcpyDeviceToHost)); 8339566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_b, d_b, sizeof(double) * N, cudaMemcpyDeviceToHost)); 8349566063dSJacob Faibussowitsch PetscCallCUDA(cudaMemcpy(h_c, d_c, sizeof(double) * N, cudaMemcpyDeviceToHost)); 835caccb7e3SMatthew G Knepley errorSTREAMkernel = STREAM_Triad_verify_double(h_b, h_c, h_a, scalar, N); 836caccb7e3SMatthew G Knepley if (errorSTREAMkernel) { 8379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " device STREAM_Triad:\t\tError detected in device STREAM_Triad, exiting\n")); 838caccb7e3SMatthew G Knepley exit(-6000); 839caccb7e3SMatthew G Knepley } 840caccb7e3SMatthew G Knepley 84119816777SMark free(h_a); 84219816777SMark free(h_b); 84319816777SMark free(h_c); 84419816777SMark } 845caccb7e3SMatthew G Knepley /* continue from here */ 84619816777SMark printResultsReadable(times,sizeof(double)); 847caccb7e3SMatthew G Knepley 848caccb7e3SMatthew G Knepley /* Free memory on device */ 8499566063dSJacob Faibussowitsch PetscCallCUDA(cudaFree(d_a)); 8509566063dSJacob Faibussowitsch PetscCallCUDA(cudaFree(d_b)); 8519566063dSJacob Faibussowitsch PetscCallCUDA(cudaFree(d_c)); 85219816777SMark 853caccb7e3SMatthew G Knepley PetscFunctionReturn(0); 854caccb7e3SMatthew G Knepley } 855caccb7e3SMatthew G Knepley 856403adfb6SMatthew G Knepley /////////////////////////////////////////////////////////////////////////// 857403adfb6SMatthew G Knepley //Print Results to Screen and File 858403adfb6SMatthew G Knepley /////////////////////////////////////////////////////////////////////////// 85919816777SMark PetscErrorCode printResultsReadable(float times[][NTIMES], const size_t bsize) 860a6dfd86eSKarl Rupp { 861403adfb6SMatthew G Knepley PetscErrorCode ierr; 862403adfb6SMatthew G Knepley PetscInt j, k; 863caccb7e3SMatthew G Knepley float avgtime[8] = {0., 0., 0., 0., 0., 0., 0., 0.}; 864caccb7e3SMatthew G Knepley float maxtime[8] = {0., 0., 0., 0., 0., 0., 0., 0.}; 865caccb7e3SMatthew G Knepley float mintime[8] = {1e30,1e30,1e30,1e30,1e30,1e30,1e30,1e30}; 86619816777SMark // char *label[8] = {"Copy: ", "Copy Opt.: ", "Scale: ", "Scale Opt: ", "Add: ", "Add Opt: ", "Triad: ", "Triad Opt: "}; 86719816777SMark const float bytes_per_kernel[8] = { 86819816777SMark 2. * bsize * N, 86919816777SMark 2. * bsize * N, 87019816777SMark 2. * bsize * N, 87119816777SMark 2. * bsize * N, 87219816777SMark 3. * bsize * N, 87319816777SMark 3. * bsize * N, 87419816777SMark 3. * bsize * N, 87519816777SMark 3. * bsize * N 876403adfb6SMatthew G Knepley }; 87719816777SMark double rate,irate; 87819816777SMark int rank,size; 879403adfb6SMatthew G Knepley PetscFunctionBegin; 8809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD,&rank)); 8819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&size)); 882403adfb6SMatthew G Knepley /* --- SUMMARY --- */ 88319816777SMark for (k = 0; k < NTIMES; ++k) { 884caccb7e3SMatthew G Knepley for (j = 0; j < 8; ++j) { 88519816777SMark avgtime[j] = avgtime[j] + (1.e-03f * times[j][k]); // millisec --> sec 886403adfb6SMatthew G Knepley mintime[j] = MIN(mintime[j], (1.e-03f * times[j][k])); 887403adfb6SMatthew G Knepley maxtime[j] = MAX(maxtime[j], (1.e-03f * times[j][k])); 888403adfb6SMatthew G Knepley } 88919816777SMark } 890caccb7e3SMatthew G Knepley for (j = 0; j < 8; ++j) { 891403adfb6SMatthew G Knepley avgtime[j] = avgtime[j]/(float)(NTIMES-1); 892403adfb6SMatthew G Knepley } 89319816777SMark j = 7; 89419816777SMark irate = 1.0E-06 * bytes_per_kernel[j]/mintime[j]; 89519816777SMark ierr = MPI_Reduce(&irate,&rate,1,MPI_DOUBLE,MPI_SUM,0,MPI_COMM_WORLD); 896dd400576SPatrick Sanan if (rank == 0) { 89719816777SMark FILE *fd; 89819816777SMark if (size == 1) { 89919816777SMark printf("%d %11.4f Rate (MB/s)\n",size, rate); 90019816777SMark fd = fopen("flops","w"); 90119816777SMark fprintf(fd,"%g\n",rate); 90219816777SMark fclose(fd); 90319816777SMark } else { 90419816777SMark double prate; 90519816777SMark fd = fopen("flops","r"); 90619816777SMark fscanf(fd,"%lg",&prate); 90719816777SMark fclose(fd); 90819816777SMark printf("%d %11.4f Rate (MB/s) %g \n", size, rate, rate/prate); 90919816777SMark } 91019816777SMark } 91119816777SMark 912403adfb6SMatthew G Knepley PetscFunctionReturn(0); 913403adfb6SMatthew G Knepley } 914