xref: /petsc/src/benchmarks/streams/CUDAVersion.cu (revision 8563dfcc547dae59d5e1dc6b8e602db41cb72b14)
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.
10403adfb6SMatthew G Knepley   The implementation is in single precision.
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 */
22403adfb6SMatthew G Knepley static char *help = "Single-Precision STREAM Benchmark implementation in CUDA\n"
23403adfb6SMatthew G Knepley                     "Performs Copy, Scale, Add, and Triad single-precision kernels\n\n";
24403adfb6SMatthew G Knepley 
25403adfb6SMatthew G Knepley #include <petscconf.h>
26403adfb6SMatthew G Knepley #include <petscsys.h>
27403adfb6SMatthew G Knepley #include <petsctime.h>
28403adfb6SMatthew G Knepley 
29403adfb6SMatthew G Knepley #define N        2000000
30caccb7e3SMatthew G Knepley #define N_DOUBLE 8000000
31403adfb6SMatthew G Knepley #define NTIMES   10
32403adfb6SMatthew G Knepley 
33403adfb6SMatthew G Knepley # ifndef MIN
34403adfb6SMatthew G Knepley # define MIN(x,y) ((x)<(y) ? (x) : (y))
35403adfb6SMatthew G Knepley # endif
36403adfb6SMatthew G Knepley # ifndef MAX
37403adfb6SMatthew G Knepley # define MAX(x,y) ((x)>(y) ? (x) : (y))
38403adfb6SMatthew G Knepley # endif
39403adfb6SMatthew G Knepley 
40403adfb6SMatthew G Knepley const float  flt_eps = 1.192092896e-07f;
41caccb7e3SMatthew G Knepley const double dbl_eps = 2.2204460492503131e-16;
42403adfb6SMatthew G Knepley 
43403adfb6SMatthew G Knepley __global__ void set_array(float *a,  float value, size_t len)
44403adfb6SMatthew G Knepley {
45403adfb6SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
46403adfb6SMatthew G Knepley   while (idx < len) {
47403adfb6SMatthew G Knepley     a[idx] = value;
48403adfb6SMatthew G Knepley     idx   += blockDim.x * gridDim.x;
49403adfb6SMatthew G Knepley   }
50403adfb6SMatthew G Knepley }
51403adfb6SMatthew G Knepley 
52caccb7e3SMatthew G Knepley __global__ void set_array_double(double *a,  double value, size_t len)
53caccb7e3SMatthew G Knepley {
54caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
55caccb7e3SMatthew G Knepley   while (idx < len) {
56caccb7e3SMatthew G Knepley     a[idx] = value;
57caccb7e3SMatthew G Knepley     idx   += blockDim.x * gridDim.x;
58caccb7e3SMatthew G Knepley   }
59caccb7e3SMatthew G Knepley }
60caccb7e3SMatthew G Knepley 
61403adfb6SMatthew G Knepley __global__ void STREAM_Copy(float *a, float *b, size_t len)
62403adfb6SMatthew G Knepley {
63403adfb6SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
64403adfb6SMatthew G Knepley   while (idx < len) {
65403adfb6SMatthew G Knepley     b[idx] = a[idx];
66403adfb6SMatthew G Knepley     idx   += blockDim.x * gridDim.x;
67403adfb6SMatthew G Knepley   }
68403adfb6SMatthew G Knepley }
69403adfb6SMatthew G Knepley 
70caccb7e3SMatthew G Knepley __global__ void STREAM_Copy_double(double *a, double *b, size_t len)
71caccb7e3SMatthew G Knepley {
72caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
73caccb7e3SMatthew G Knepley   while (idx < len) {
74caccb7e3SMatthew G Knepley     b[idx] = a[idx];
75caccb7e3SMatthew G Knepley     idx   += blockDim.x * gridDim.x;
76caccb7e3SMatthew G Knepley   }
77caccb7e3SMatthew G Knepley }
78caccb7e3SMatthew G Knepley 
79403adfb6SMatthew G Knepley __global__ void STREAM_Copy_Optimized(float *a, float *b, size_t len)
80403adfb6SMatthew G Knepley {
81403adfb6SMatthew G Knepley   /*
82403adfb6SMatthew G Knepley    * Ensure size of thread index space is as large as or greater than
83403adfb6SMatthew G Knepley    * vector index space else return.
84403adfb6SMatthew G Knepley    */
85403adfb6SMatthew G Knepley   if (blockDim.x * gridDim.x < len) return;
86403adfb6SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
87403adfb6SMatthew G Knepley   if (idx < len) b[idx] = a[idx];
88403adfb6SMatthew G Knepley }
89403adfb6SMatthew G Knepley 
90caccb7e3SMatthew G Knepley __global__ void STREAM_Copy_Optimized_double(double *a, double *b, size_t len)
91caccb7e3SMatthew G Knepley {
92caccb7e3SMatthew G Knepley   /*
93caccb7e3SMatthew G Knepley    * Ensure size of thread index space is as large as or greater than
94caccb7e3SMatthew G Knepley    * vector index space else return.
95caccb7e3SMatthew G Knepley    */
96caccb7e3SMatthew G Knepley   if (blockDim.x * gridDim.x < len) return;
97caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
98caccb7e3SMatthew G Knepley   if (idx < len) b[idx] = a[idx];
99caccb7e3SMatthew G Knepley }
100caccb7e3SMatthew G Knepley 
101403adfb6SMatthew G Knepley __global__ void STREAM_Scale(float *a, float *b, float scale,  size_t len)
102403adfb6SMatthew G Knepley {
103403adfb6SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
104403adfb6SMatthew G Knepley   while (idx < len) {
105403adfb6SMatthew G Knepley     b[idx] = scale* a[idx];
106403adfb6SMatthew G Knepley     idx   += blockDim.x * gridDim.x;
107403adfb6SMatthew G Knepley   }
108403adfb6SMatthew G Knepley }
109403adfb6SMatthew G Knepley 
110caccb7e3SMatthew G Knepley __global__ void STREAM_Scale_double(double *a, double *b, double scale,  size_t len)
111caccb7e3SMatthew G Knepley {
112caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
113caccb7e3SMatthew G Knepley   while (idx < len) {
114caccb7e3SMatthew G Knepley     b[idx] = scale* a[idx];
115caccb7e3SMatthew G Knepley     idx   += blockDim.x * gridDim.x;
116caccb7e3SMatthew G Knepley   }
117caccb7e3SMatthew G Knepley }
118caccb7e3SMatthew G Knepley 
119caccb7e3SMatthew G Knepley __global__ void STREAM_Scale_Optimized(float *a, float *b, float scale,  size_t len)
120caccb7e3SMatthew G Knepley {
121caccb7e3SMatthew G Knepley   /*
122caccb7e3SMatthew G Knepley    * Ensure size of thread index space is as large as or greater than
123caccb7e3SMatthew G Knepley    * vector index space else return.
124caccb7e3SMatthew G Knepley    */
125caccb7e3SMatthew G Knepley   if (blockDim.x * gridDim.x < len) return;
126caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
127caccb7e3SMatthew G Knepley   if (idx < len) b[idx] = scale* a[idx];
128caccb7e3SMatthew G Knepley }
129caccb7e3SMatthew G Knepley 
130caccb7e3SMatthew G Knepley __global__ void STREAM_Scale_Optimized_double(double *a, double *b, double scale,  size_t len)
131caccb7e3SMatthew G Knepley {
132caccb7e3SMatthew G Knepley   /*
133caccb7e3SMatthew G Knepley    * Ensure size of thread index space is as large as or greater than
134caccb7e3SMatthew G Knepley    * vector index space else return.
135caccb7e3SMatthew G Knepley    */
136caccb7e3SMatthew G Knepley   if (blockDim.x * gridDim.x < len) return;
137caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
138caccb7e3SMatthew G Knepley   if (idx < len) b[idx] = scale* a[idx];
139caccb7e3SMatthew G Knepley }
140caccb7e3SMatthew G Knepley 
141403adfb6SMatthew G Knepley __global__ void STREAM_Add(float *a, float *b, float *c,  size_t len)
142403adfb6SMatthew G Knepley {
143403adfb6SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
144403adfb6SMatthew G Knepley   while (idx < len) {
145403adfb6SMatthew G Knepley     c[idx] = a[idx]+b[idx];
146403adfb6SMatthew G Knepley     idx   += blockDim.x * gridDim.x;
147403adfb6SMatthew G Knepley   }
148403adfb6SMatthew G Knepley }
149403adfb6SMatthew G Knepley 
150caccb7e3SMatthew G Knepley __global__ void STREAM_Add_double(double *a, double *b, double *c,  size_t len)
151caccb7e3SMatthew G Knepley {
152caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
153caccb7e3SMatthew G Knepley   while (idx < len) {
154caccb7e3SMatthew G Knepley     c[idx] = a[idx]+b[idx];
155caccb7e3SMatthew G Knepley     idx   += blockDim.x * gridDim.x;
156caccb7e3SMatthew G Knepley   }
157caccb7e3SMatthew G Knepley }
158caccb7e3SMatthew G Knepley 
159caccb7e3SMatthew G Knepley __global__ void STREAM_Add_Optimized(float *a, float *b, float *c,  size_t len)
160caccb7e3SMatthew G Knepley {
161caccb7e3SMatthew G Knepley   /*
162caccb7e3SMatthew G Knepley    * Ensure size of thread index space is as large as or greater than
163caccb7e3SMatthew G Knepley    * vector index space else return.
164caccb7e3SMatthew G Knepley    */
165caccb7e3SMatthew G Knepley   if (blockDim.x * gridDim.x < len) return;
166caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
167caccb7e3SMatthew G Knepley   if (idx < len) c[idx] = a[idx]+b[idx];
168caccb7e3SMatthew G Knepley }
169caccb7e3SMatthew G Knepley 
170caccb7e3SMatthew G Knepley __global__ void STREAM_Add_Optimized_double(double *a, double *b, double *c,  size_t len)
171caccb7e3SMatthew G Knepley {
172caccb7e3SMatthew G Knepley   /*
173caccb7e3SMatthew G Knepley    * Ensure size of thread index space is as large as or greater than
174caccb7e3SMatthew G Knepley    * vector index space else return.
175caccb7e3SMatthew G Knepley    */
176caccb7e3SMatthew G Knepley   if (blockDim.x * gridDim.x < len) return;
177caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
178caccb7e3SMatthew G Knepley   if (idx < len) c[idx] = a[idx]+b[idx];
179caccb7e3SMatthew G Knepley }
180caccb7e3SMatthew G Knepley 
181403adfb6SMatthew G Knepley __global__ void STREAM_Triad(float *a, float *b, float *c, float scalar, size_t len)
182403adfb6SMatthew G Knepley {
183403adfb6SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
184403adfb6SMatthew G Knepley   while (idx < len) {
185403adfb6SMatthew G Knepley     c[idx] = a[idx]+scalar*b[idx];
186403adfb6SMatthew G Knepley     idx   += blockDim.x * gridDim.x;
187403adfb6SMatthew G Knepley   }
188403adfb6SMatthew G Knepley }
189403adfb6SMatthew G Knepley 
190caccb7e3SMatthew G Knepley __global__ void STREAM_Triad_double(double *a, double *b, double *c, double scalar, size_t len)
191caccb7e3SMatthew G Knepley {
192caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
193caccb7e3SMatthew G Knepley   while (idx < len) {
194caccb7e3SMatthew G Knepley     c[idx] = a[idx]+scalar*b[idx];
195caccb7e3SMatthew G Knepley     idx   += blockDim.x * gridDim.x;
196caccb7e3SMatthew G Knepley   }
197caccb7e3SMatthew G Knepley }
198caccb7e3SMatthew G Knepley 
199caccb7e3SMatthew G Knepley __global__ void STREAM_Triad_Optimized(float *a, float *b, float *c, float scalar, size_t len)
200caccb7e3SMatthew G Knepley {
201caccb7e3SMatthew G Knepley   /*
202caccb7e3SMatthew G Knepley    * Ensure size of thread index space is as large as or greater than
203caccb7e3SMatthew G Knepley    * vector index space else return.
204caccb7e3SMatthew G Knepley    */
205caccb7e3SMatthew G Knepley   if (blockDim.x * gridDim.x < len) return;
206caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
207caccb7e3SMatthew G Knepley   if (idx < len) c[idx] = a[idx]+scalar*b[idx];
208caccb7e3SMatthew G Knepley }
209caccb7e3SMatthew G Knepley 
210caccb7e3SMatthew G Knepley __global__ void STREAM_Triad_Optimized_double(double *a, double *b, double *c, double scalar, size_t len)
211caccb7e3SMatthew G Knepley {
212caccb7e3SMatthew G Knepley   /*
213caccb7e3SMatthew G Knepley    * Ensure size of thread index space is as large as or greater than
214caccb7e3SMatthew G Knepley    * vector index space else return.
215caccb7e3SMatthew G Knepley    */
216caccb7e3SMatthew G Knepley   if (blockDim.x * gridDim.x < len) return;
217caccb7e3SMatthew G Knepley   size_t idx = threadIdx.x + blockIdx.x * blockDim.x;
218caccb7e3SMatthew G Knepley   if (idx < len) c[idx] = a[idx]+scalar*b[idx];
219caccb7e3SMatthew G Knepley }
220caccb7e3SMatthew G Knepley 
221403adfb6SMatthew G Knepley /* Host side verification routines */
222a6dfd86eSKarl Rupp bool STREAM_Copy_verify(float *a, float *b, size_t len)
223a6dfd86eSKarl Rupp {
224403adfb6SMatthew G Knepley   size_t idx;
225403adfb6SMatthew G Knepley   bool   bDifferent = false;
226403adfb6SMatthew G Knepley 
227403adfb6SMatthew G Knepley   for (idx = 0; idx < len && !bDifferent; idx++) {
228403adfb6SMatthew G Knepley     float expectedResult     = a[idx];
229403adfb6SMatthew G Knepley     float diffResultExpected = (b[idx] - expectedResult);
230403adfb6SMatthew G Knepley     float relErrorULPS       = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps;
231403adfb6SMatthew G Knepley     /* element-wise relative error determination */
232403adfb6SMatthew G Knepley     bDifferent = (relErrorULPS > 2.f);
233403adfb6SMatthew G Knepley   }
234403adfb6SMatthew G Knepley 
235403adfb6SMatthew G Knepley   return bDifferent;
236403adfb6SMatthew G Knepley }
237403adfb6SMatthew G Knepley 
238a6dfd86eSKarl Rupp bool STREAM_Copy_verify_double(double *a, double *b, size_t len)
239a6dfd86eSKarl Rupp {
240caccb7e3SMatthew G Knepley   size_t idx;
241caccb7e3SMatthew G Knepley   bool   bDifferent = false;
242caccb7e3SMatthew G Knepley 
243caccb7e3SMatthew G Knepley   for (idx = 0; idx < len && !bDifferent; idx++) {
244caccb7e3SMatthew G Knepley     double expectedResult     = a[idx];
245caccb7e3SMatthew G Knepley     double diffResultExpected = (b[idx] - expectedResult);
246caccb7e3SMatthew G Knepley     double relErrorULPS       = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps;
247caccb7e3SMatthew G Knepley     /* element-wise relative error determination */
248caccb7e3SMatthew G Knepley     bDifferent = (relErrorULPS > 2.);
249caccb7e3SMatthew G Knepley   }
250caccb7e3SMatthew G Knepley 
251caccb7e3SMatthew G Knepley   return bDifferent;
252caccb7e3SMatthew G Knepley }
253caccb7e3SMatthew G Knepley 
254a6dfd86eSKarl Rupp bool STREAM_Scale_verify(float *a, float *b, float scale, size_t len)
255a6dfd86eSKarl Rupp {
256403adfb6SMatthew G Knepley   size_t idx;
257403adfb6SMatthew G Knepley   bool   bDifferent = false;
258403adfb6SMatthew G Knepley 
259403adfb6SMatthew G Knepley   for (idx = 0; idx < len && !bDifferent; idx++) {
260403adfb6SMatthew G Knepley     float expectedResult     = scale*a[idx];
261403adfb6SMatthew G Knepley     float diffResultExpected = (b[idx] - expectedResult);
262403adfb6SMatthew G Knepley     float relErrorULPS       = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps;
263403adfb6SMatthew G Knepley     /* element-wise relative error determination */
264403adfb6SMatthew G Knepley     bDifferent = (relErrorULPS > 2.f);
265403adfb6SMatthew G Knepley   }
266403adfb6SMatthew G Knepley 
267403adfb6SMatthew G Knepley   return bDifferent;
268403adfb6SMatthew G Knepley }
269403adfb6SMatthew G Knepley 
270a6dfd86eSKarl Rupp bool STREAM_Scale_verify_double(double *a, double *b, double scale, size_t len)
271a6dfd86eSKarl Rupp {
272caccb7e3SMatthew G Knepley   size_t idx;
273caccb7e3SMatthew G Knepley   bool   bDifferent = false;
274caccb7e3SMatthew G Knepley 
275caccb7e3SMatthew G Knepley   for (idx = 0; idx < len && !bDifferent; idx++) {
276caccb7e3SMatthew G Knepley     double expectedResult     = scale*a[idx];
277caccb7e3SMatthew G Knepley     double diffResultExpected = (b[idx] - expectedResult);
278caccb7e3SMatthew G Knepley     double relErrorULPS       = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps;
279caccb7e3SMatthew G Knepley     /* element-wise relative error determination */
280caccb7e3SMatthew G Knepley     bDifferent = (relErrorULPS > 2.);
281caccb7e3SMatthew G Knepley   }
282caccb7e3SMatthew G Knepley 
283caccb7e3SMatthew G Knepley   return bDifferent;
284caccb7e3SMatthew G Knepley }
285caccb7e3SMatthew G Knepley 
286a6dfd86eSKarl Rupp bool STREAM_Add_verify(float *a, float *b, float *c, size_t len)
287a6dfd86eSKarl Rupp {
288403adfb6SMatthew G Knepley   size_t idx;
289403adfb6SMatthew G Knepley   bool   bDifferent = false;
290403adfb6SMatthew G Knepley 
291403adfb6SMatthew G Knepley   for (idx = 0; idx < len && !bDifferent; idx++) {
292403adfb6SMatthew G Knepley     float expectedResult     = a[idx] + b[idx];
293403adfb6SMatthew G Knepley     float diffResultExpected = (c[idx] - expectedResult);
294403adfb6SMatthew G Knepley     float relErrorULPS       = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps;
295403adfb6SMatthew G Knepley     /* element-wise relative error determination */
296403adfb6SMatthew G Knepley     bDifferent = (relErrorULPS > 2.f);
297403adfb6SMatthew G Knepley   }
298403adfb6SMatthew G Knepley 
299403adfb6SMatthew G Knepley   return bDifferent;
300403adfb6SMatthew G Knepley }
301403adfb6SMatthew G Knepley 
302a6dfd86eSKarl Rupp bool STREAM_Add_verify_double(double *a, double *b, double *c, size_t len)
303a6dfd86eSKarl Rupp {
304caccb7e3SMatthew G Knepley   size_t idx;
305caccb7e3SMatthew G Knepley   bool   bDifferent = false;
306caccb7e3SMatthew G Knepley 
307caccb7e3SMatthew G Knepley   for (idx = 0; idx < len && !bDifferent; idx++) {
308caccb7e3SMatthew G Knepley     double expectedResult     = a[idx] + b[idx];
309caccb7e3SMatthew G Knepley     double diffResultExpected = (c[idx] - expectedResult);
310caccb7e3SMatthew G Knepley     double relErrorULPS       = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps;
311caccb7e3SMatthew G Knepley     /* element-wise relative error determination */
312caccb7e3SMatthew G Knepley     bDifferent = (relErrorULPS > 2.);
313caccb7e3SMatthew G Knepley   }
314caccb7e3SMatthew G Knepley 
315caccb7e3SMatthew G Knepley   return bDifferent;
316caccb7e3SMatthew G Knepley }
317caccb7e3SMatthew G Knepley 
318a6dfd86eSKarl Rupp bool STREAM_Triad_verify(float *a, float *b, float *c, float scalar, size_t len)
319a6dfd86eSKarl Rupp {
320403adfb6SMatthew G Knepley   size_t idx;
321403adfb6SMatthew G Knepley   bool   bDifferent = false;
322403adfb6SMatthew G Knepley 
323403adfb6SMatthew G Knepley   for (idx = 0; idx < len && !bDifferent; idx++) {
324403adfb6SMatthew G Knepley     float expectedResult     = a[idx] + scalar*b[idx];
325403adfb6SMatthew G Knepley     float diffResultExpected = (c[idx] - expectedResult);
326403adfb6SMatthew G Knepley     float relErrorULPS       = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps;
327403adfb6SMatthew G Knepley     /* element-wise relative error determination */
328403adfb6SMatthew G Knepley     bDifferent = (relErrorULPS > 3.f);
329403adfb6SMatthew G Knepley   }
330403adfb6SMatthew G Knepley 
331403adfb6SMatthew G Knepley   return bDifferent;
332403adfb6SMatthew G Knepley }
333403adfb6SMatthew G Knepley 
334a6dfd86eSKarl Rupp bool STREAM_Triad_verify_double(double *a, double *b, double *c, double scalar, size_t len)
335a6dfd86eSKarl Rupp {
336caccb7e3SMatthew G Knepley   size_t idx;
337caccb7e3SMatthew G Knepley   bool   bDifferent = false;
338caccb7e3SMatthew G Knepley 
339caccb7e3SMatthew G Knepley   for (idx = 0; idx < len && !bDifferent; idx++) {
340caccb7e3SMatthew G Knepley     double expectedResult     = a[idx] + scalar*b[idx];
341caccb7e3SMatthew G Knepley     double diffResultExpected = (c[idx] - expectedResult);
342caccb7e3SMatthew G Knepley     double relErrorULPS       = (fabsf(diffResultExpected)/fabsf(expectedResult))/flt_eps;
343caccb7e3SMatthew G Knepley     /* element-wise relative error determination */
344caccb7e3SMatthew G Knepley     bDifferent = (relErrorULPS > 3.);
345caccb7e3SMatthew G Knepley   }
346caccb7e3SMatthew G Knepley 
347caccb7e3SMatthew G Knepley   return bDifferent;
348caccb7e3SMatthew G Knepley }
349caccb7e3SMatthew G Knepley 
350403adfb6SMatthew G Knepley /* forward declarations */
351caccb7e3SMatthew G Knepley PetscErrorCode setupStream(PetscInt device, PetscBool runDouble, PetscBool cpuTiming);
352403adfb6SMatthew G Knepley PetscErrorCode runStream(const PetscInt iNumThreadsPerBlock, PetscBool bDontUseGPUTiming);
353caccb7e3SMatthew G Knepley PetscErrorCode runStreamDouble(const PetscInt iNumThreadsPerBlock, PetscBool bDontUseGPUTiming);
354403adfb6SMatthew G Knepley PetscErrorCode printResultsReadable(float times[][NTIMES]);
355403adfb6SMatthew G Knepley 
356403adfb6SMatthew G Knepley int main(int argc, char *argv[])
357403adfb6SMatthew G Knepley {
358403adfb6SMatthew G Knepley   PetscInt       device    = 0;
359caccb7e3SMatthew G Knepley   PetscBool      runDouble = PETSC_FALSE;
360403adfb6SMatthew G Knepley   PetscBool      cpuTiming = PETSC_FALSE;
361403adfb6SMatthew G Knepley   PetscErrorCode ierr;
362403adfb6SMatthew G Knepley 
363403adfb6SMatthew G Knepley   ierr = PetscInitialize(&argc, &argv, 0, help);CHKERRQ(ierr);
364caccb7e3SMatthew G Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, "[Single and Double-Precision Device-Only STREAM Benchmark implementation in CUDA]\n");CHKERRQ(ierr);
365403adfb6SMatthew G Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, "%s Starting...\n\n", argv[0]);CHKERRQ(ierr);
366403adfb6SMatthew G Knepley 
367403adfb6SMatthew G Knepley   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "STREAM Benchmark Options", "STREAM");CHKERRQ(ierr);
3680298fd71SBarry Smith   ierr = PetscOptionsInt("-device", "Specify the CUDA device to be used", "STREAM", device, &device, NULL);CHKERRQ(ierr);
3690298fd71SBarry Smith   ierr = PetscOptionsBool("-double",    "Also run double precision tests",   "STREAM", runDouble, &runDouble, NULL);CHKERRQ(ierr);
3700298fd71SBarry Smith   ierr = PetscOptionsBool("-cputiming", "Force CPU-based timing to be used", "STREAM", cpuTiming, &cpuTiming, NULL);CHKERRQ(ierr);
371403adfb6SMatthew G Knepley   ierr = PetscOptionsEnd();
372403adfb6SMatthew G Knepley 
373caccb7e3SMatthew G Knepley   ierr = setupStream(device, runDouble, cpuTiming);
374403adfb6SMatthew G Knepley   if (ierr >= 0) {
375403adfb6SMatthew G Knepley     PetscErrorCode ierr2 = PetscPrintf(PETSC_COMM_SELF, "\n[streamBenchmark] - results:\t%s\n\n", (ierr == 0) ? "PASSES" : "FAILED");CHKERRQ(ierr2);
376403adfb6SMatthew G Knepley   }
377403adfb6SMatthew G Knepley   PetscFinalize();
378403adfb6SMatthew G Knepley   return 0;
379403adfb6SMatthew G Knepley }
380403adfb6SMatthew G Knepley 
381403adfb6SMatthew G Knepley ///////////////////////////////////////////////////////////////////////////////
382403adfb6SMatthew G Knepley //Run the appropriate tests
383403adfb6SMatthew G Knepley ///////////////////////////////////////////////////////////////////////////////
384caccb7e3SMatthew G Knepley PetscErrorCode setupStream(PetscInt deviceNum, PetscBool runDouble, PetscBool cpuTiming)
385403adfb6SMatthew G Knepley {
386403adfb6SMatthew G Knepley   PetscInt       iNumThreadsPerBlock = 128;
387403adfb6SMatthew G Knepley   PetscErrorCode ierr;
388403adfb6SMatthew G Knepley 
389403adfb6SMatthew G Knepley   PetscFunctionBegin;
390403adfb6SMatthew G Knepley   // Check device
391403adfb6SMatthew G Knepley   {
392403adfb6SMatthew G Knepley     int deviceCount;
393403adfb6SMatthew G Knepley 
394403adfb6SMatthew G Knepley     cudaGetDeviceCount(&deviceCount);
395403adfb6SMatthew G Knepley     if (deviceCount == 0) {
396403adfb6SMatthew G Knepley       ierr = PetscPrintf(PETSC_COMM_SELF, "!!!!!No devices found!!!!!\n");CHKERRQ(ierr);
397403adfb6SMatthew G Knepley       return -1000;
398403adfb6SMatthew G Knepley     }
399403adfb6SMatthew G Knepley 
400403adfb6SMatthew G Knepley     if (deviceNum >= deviceCount || deviceNum < 0) {
401403adfb6SMatthew G Knepley       ierr      = PetscPrintf(PETSC_COMM_SELF, "\n!!!!!Invalid GPU number %d given hence default gpu %d will be used !!!!!\n", deviceNum, 0);CHKERRQ(ierr);
402403adfb6SMatthew G Knepley       deviceNum = 0;
403403adfb6SMatthew G Knepley     }
404403adfb6SMatthew G Knepley   }
405403adfb6SMatthew G Knepley 
406403adfb6SMatthew G Knepley   cudaSetDevice(deviceNum);
407403adfb6SMatthew G Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, "Running on...\n\n");CHKERRQ(ierr);
408403adfb6SMatthew G Knepley   cudaDeviceProp deviceProp;
409403adfb6SMatthew G Knepley   if (cudaGetDeviceProperties(&deviceProp, deviceNum) == cudaSuccess) {
410403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " Device %d: %s\n", deviceNum, deviceProp.name);CHKERRQ(ierr);
411403adfb6SMatthew G Knepley   } else {
412403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " Unable to determine device %d properties, exiting\n");CHKERRQ(ierr);
413403adfb6SMatthew G Knepley     return -1;
414403adfb6SMatthew G Knepley   }
415403adfb6SMatthew G Knepley 
416caccb7e3SMatthew G Knepley   if (runDouble && deviceProp.major == 1 && deviceProp.minor < 3) {
417caccb7e3SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " Unable to run double-precision STREAM benchmark on a compute capability GPU less than 1.3\n");CHKERRQ(ierr);
418caccb7e3SMatthew G Knepley     return -1;
419caccb7e3SMatthew G Knepley   }
4206f2b61bcSKarl Rupp   if (deviceProp.major == 2 && deviceProp.minor == 1) iNumThreadsPerBlock = 192; /* GF104 architecture / 48 CUDA Cores per MP */
4216f2b61bcSKarl Rupp   else iNumThreadsPerBlock = 128; /* GF100 architecture / 32 CUDA Cores per MP */
422403adfb6SMatthew G Knepley 
423403adfb6SMatthew G Knepley   if (cpuTiming) {
424403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " Using cpu-only timer.\n");CHKERRQ(ierr);
425403adfb6SMatthew G Knepley   }
426403adfb6SMatthew G Knepley 
427403adfb6SMatthew G Knepley   ierr = runStream(iNumThreadsPerBlock, cpuTiming);CHKERRQ(ierr);
428caccb7e3SMatthew G Knepley   if (runDouble) {
429caccb7e3SMatthew G Knepley     ierr = cudaSetDeviceFlags(cudaDeviceBlockingSync);CHKERRQ(ierr);
430caccb7e3SMatthew G Knepley     ierr = runStreamDouble(iNumThreadsPerBlock, cpuTiming);CHKERRQ(ierr);
431caccb7e3SMatthew G Knepley   }
432403adfb6SMatthew G Knepley   PetscFunctionReturn(0);
433403adfb6SMatthew G Knepley }
434403adfb6SMatthew G Knepley 
435403adfb6SMatthew G Knepley ///////////////////////////////////////////////////////////////////////////
436403adfb6SMatthew G Knepley // runStream
437403adfb6SMatthew G Knepley ///////////////////////////////////////////////////////////////////////////
438403adfb6SMatthew G Knepley PetscErrorCode runStream(const PetscInt iNumThreadsPerBlock, PetscBool bDontUseGPUTiming)
439403adfb6SMatthew G Knepley {
440403adfb6SMatthew G Knepley   float          *d_a, *d_b, *d_c;
441403adfb6SMatthew G Knepley   int            k;
442caccb7e3SMatthew G Knepley   float          times[8][NTIMES];
443403adfb6SMatthew G Knepley   float          scalar;
444403adfb6SMatthew G Knepley   PetscErrorCode ierr;
445403adfb6SMatthew G Knepley 
446403adfb6SMatthew G Knepley   PetscFunctionBegin;
447403adfb6SMatthew G Knepley   /* Allocate memory on device */
448403adfb6SMatthew G Knepley   ierr = cudaMalloc((void**)&d_a, sizeof(float)*N);CHKERRQ(ierr);
449403adfb6SMatthew G Knepley   ierr = cudaMalloc((void**)&d_b, sizeof(float)*N);CHKERRQ(ierr);
450403adfb6SMatthew G Knepley   ierr = cudaMalloc((void**)&d_c, sizeof(float)*N);CHKERRQ(ierr);
451403adfb6SMatthew G Knepley 
452403adfb6SMatthew G Knepley   /* Compute execution configuration */
453403adfb6SMatthew G Knepley 
454403adfb6SMatthew G Knepley   dim3 dimBlock(iNumThreadsPerBlock); /* (iNumThreadsPerBlock,1,1) */
455403adfb6SMatthew G Knepley   dim3 dimGrid(N/dimBlock.x); /* (N/dimBlock.x,1,1) */
456403adfb6SMatthew G Knepley   if (N % dimBlock.x != 0) dimGrid.x+=1;
457403adfb6SMatthew G Knepley 
458403adfb6SMatthew G Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, " Array size (single precision) = %u\n",N);CHKERRQ(ierr);
459403adfb6SMatthew G Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, " using %u threads per block, %u blocks\n",dimBlock.x,dimGrid.x);CHKERRQ(ierr);
460403adfb6SMatthew G Knepley 
461403adfb6SMatthew G Knepley   /* Initialize memory on the device */
462403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_a, 2.f, N);
463403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N);
464403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N);
465403adfb6SMatthew G Knepley 
466403adfb6SMatthew G Knepley   /* --- MAIN LOOP --- repeat test cases NTIMES times --- */
467403adfb6SMatthew G Knepley   PetscLogDouble cpuTimer = 0.0;
468403adfb6SMatthew G Knepley   cudaEvent_t    start, stop;
469403adfb6SMatthew G Knepley 
470403adfb6SMatthew G Knepley   /* both timers report msec */
471403adfb6SMatthew G Knepley   ierr = cudaEventCreate(&start);CHKERRQ(ierr); /* gpu timer facility */
472403adfb6SMatthew G Knepley   ierr = cudaEventCreate(&stop);CHKERRQ(ierr);  /* gpu timer facility */
473403adfb6SMatthew G Knepley 
474403adfb6SMatthew G Knepley   scalar=3.0f;
475403adfb6SMatthew G Knepley   for (k = 0; k < NTIMES; ++k) {
476*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
477403adfb6SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
478403adfb6SMatthew G Knepley     STREAM_Copy<<<dimGrid,dimBlock>>>(d_a, d_c, N);
479403adfb6SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
480403adfb6SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
481403adfb6SMatthew G Knepley     //get the the total elapsed time in ms
482*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
4836f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[0][k] = cpuTimer;
4846f2b61bcSKarl Rupp     else {
485403adfb6SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[0][k], start, stop);CHKERRQ(ierr);
486403adfb6SMatthew G Knepley     }
487403adfb6SMatthew G Knepley 
488403adfb6SMatthew G Knepley     cpuTimer = 0.0;
489*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
490403adfb6SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
491403adfb6SMatthew G Knepley     STREAM_Copy_Optimized<<<dimGrid,dimBlock>>>(d_a, d_c, N);
492403adfb6SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
493403adfb6SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
494403adfb6SMatthew G Knepley     //get the the total elapsed time in ms
495*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
4966f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[1][k] = cpuTimer;
4976f2b61bcSKarl Rupp     else {
498403adfb6SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[1][k], start, stop);CHKERRQ(ierr);
499403adfb6SMatthew G Knepley     }
500403adfb6SMatthew G Knepley 
501403adfb6SMatthew G Knepley     cpuTimer = 0.0;
502*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
503403adfb6SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
504403adfb6SMatthew G Knepley     STREAM_Scale<<<dimGrid,dimBlock>>>(d_b, d_c, scalar,  N);
505403adfb6SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
506403adfb6SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
507403adfb6SMatthew G Knepley     //get the the total elapsed time in ms
508*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
5096f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[2][k] = cpuTimer;
5106f2b61bcSKarl Rupp     else {
511403adfb6SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[2][k], start, stop);CHKERRQ(ierr);
512403adfb6SMatthew G Knepley     }
513403adfb6SMatthew G Knepley 
514403adfb6SMatthew G Knepley     cpuTimer = 0.0;
515*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
516403adfb6SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
517caccb7e3SMatthew G Knepley     STREAM_Scale_Optimized<<<dimGrid,dimBlock>>>(d_b, d_c, scalar,  N);
518403adfb6SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
519403adfb6SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
520403adfb6SMatthew G Knepley     //get the the total elapsed time in ms
521*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
5226f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[3][k] = cpuTimer;
5236f2b61bcSKarl Rupp     else {
524403adfb6SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[3][k], start, stop);CHKERRQ(ierr);
525403adfb6SMatthew G Knepley     }
526403adfb6SMatthew G Knepley 
527403adfb6SMatthew G Knepley     cpuTimer = 0.0;
528*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
529403adfb6SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
530caccb7e3SMatthew G Knepley     STREAM_Add<<<dimGrid,dimBlock>>>(d_a, d_b, d_c,  N);
531403adfb6SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
532403adfb6SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
533403adfb6SMatthew G Knepley     //get the the total elapsed time in ms
534*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
5356f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[4][k] = cpuTimer;
5366f2b61bcSKarl Rupp     else {
537403adfb6SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[4][k], start, stop);CHKERRQ(ierr);
538403adfb6SMatthew G Knepley     }
539403adfb6SMatthew G Knepley 
540caccb7e3SMatthew G Knepley     cpuTimer = 0.0;
541*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
542caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
543caccb7e3SMatthew G Knepley     STREAM_Add_Optimized<<<dimGrid,dimBlock>>>(d_a, d_b, d_c,  N);
544caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
545caccb7e3SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
546caccb7e3SMatthew G Knepley     //get the the total elapsed time in ms
547*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
5486f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[5][k] = cpuTimer;
5496f2b61bcSKarl Rupp     else {
550caccb7e3SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[5][k], start, stop);CHKERRQ(ierr);
551caccb7e3SMatthew G Knepley     }
552caccb7e3SMatthew G Knepley 
553caccb7e3SMatthew G Knepley     cpuTimer = 0.0;
554*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
555caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
556caccb7e3SMatthew G Knepley     STREAM_Triad<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar,  N);
557caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
558caccb7e3SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
559caccb7e3SMatthew G Knepley     //get the the total elapsed time in ms
560*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
5616f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[6][k] = cpuTimer;
5626f2b61bcSKarl Rupp     else {
563caccb7e3SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[6][k], start, stop);CHKERRQ(ierr);
564caccb7e3SMatthew G Knepley     }
565caccb7e3SMatthew G Knepley 
566caccb7e3SMatthew G Knepley     cpuTimer = 0.0;
567*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
568caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
569caccb7e3SMatthew G Knepley     STREAM_Triad_Optimized<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar,  N);
570caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
571caccb7e3SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
572caccb7e3SMatthew G Knepley     //get the the total elapsed time in ms
573*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
5746f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[7][k] = cpuTimer;
5756f2b61bcSKarl Rupp     else {
576caccb7e3SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[7][k], start, stop);CHKERRQ(ierr);
577caccb7e3SMatthew G Knepley     }
578caccb7e3SMatthew G Knepley 
579403adfb6SMatthew G Knepley   }
580403adfb6SMatthew G Knepley 
581403adfb6SMatthew G Knepley   /* verify kernels */
582403adfb6SMatthew G Knepley   float *h_a, *h_b, *h_c;
583403adfb6SMatthew G Knepley   bool  errorSTREAMkernel = true;
584403adfb6SMatthew G Knepley 
585403adfb6SMatthew G Knepley   if ((h_a = (float*)calloc(N, sizeof(float))) == (float*)NULL) {
586403adfb6SMatthew G Knepley     printf("Unable to allocate array h_a, exiting ...\n");
587403adfb6SMatthew G Knepley     exit(1);
588403adfb6SMatthew G Knepley   }
589403adfb6SMatthew G Knepley   if ((h_b = (float*)calloc(N, sizeof(float))) == (float*)NULL) {
590403adfb6SMatthew G Knepley     printf("Unable to allocate array h_b, exiting ...\n");
591403adfb6SMatthew G Knepley     exit(1);
592403adfb6SMatthew G Knepley   }
593403adfb6SMatthew G Knepley 
594403adfb6SMatthew G Knepley   if ((h_c = (float*)calloc(N, sizeof(float))) == (float*)NULL) {
595403adfb6SMatthew G Knepley     printf("Unalbe to allocate array h_c, exiting ...\n");
596403adfb6SMatthew G Knepley     exit(1);
597403adfb6SMatthew G Knepley   }
598403adfb6SMatthew G Knepley 
599403adfb6SMatthew G Knepley   /*
600403adfb6SMatthew G Knepley    * perform kernel, copy device memory into host memory and verify each
601403adfb6SMatthew G Knepley    * device kernel output
602403adfb6SMatthew G Knepley    */
603403adfb6SMatthew G Knepley 
604403adfb6SMatthew G Knepley   /* Initialize memory on the device */
605403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_a, 2.f, N);
606403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N);
607403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N);
608403adfb6SMatthew G Knepley 
609403adfb6SMatthew G Knepley   STREAM_Copy<<<dimGrid,dimBlock>>>(d_a, d_c, N);
610403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_a, d_a, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
611403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_c, d_c, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
612403adfb6SMatthew G Knepley   errorSTREAMkernel = STREAM_Copy_verify(h_a, h_c, N);
613403adfb6SMatthew G Knepley   if (errorSTREAMkernel) {
614403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy:\t\tError detected in device STREAM_Copy, exiting\n");CHKERRQ(ierr);
615403adfb6SMatthew G Knepley     exit(-2000);
616403adfb6SMatthew G Knepley   } else {
617403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy:\t\tPass\n");CHKERRQ(ierr);
618403adfb6SMatthew G Knepley   }
619403adfb6SMatthew G Knepley 
620403adfb6SMatthew G Knepley   /* Initialize memory on the device */
621403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_a, 2.f, N);
622403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N);
623403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N);
624403adfb6SMatthew G Knepley 
625403adfb6SMatthew G Knepley   STREAM_Copy_Optimized<<<dimGrid,dimBlock>>>(d_a, d_c, N);
626403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_a, d_a, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
627403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_c, d_c, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
628403adfb6SMatthew G Knepley   errorSTREAMkernel = STREAM_Copy_verify(h_a, h_c, N);
629403adfb6SMatthew G Knepley   if (errorSTREAMkernel) {
630403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy_Optimized:\tError detected in device STREAM_Copy_Optimized, exiting\n");CHKERRQ(ierr);
631403adfb6SMatthew G Knepley     exit(-3000);
632403adfb6SMatthew G Knepley   } else {
633403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy_Optimized:\tPass\n");CHKERRQ(ierr);
634403adfb6SMatthew G Knepley   }
635403adfb6SMatthew G Knepley 
636403adfb6SMatthew G Knepley   /* Initialize memory on the device */
637403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N);
638403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N);
639403adfb6SMatthew G Knepley 
640403adfb6SMatthew G Knepley   STREAM_Scale<<<dimGrid,dimBlock>>>(d_b, d_c, scalar, N);
641403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_b, d_b, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
642403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_c, d_c, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
643403adfb6SMatthew G Knepley   errorSTREAMkernel = STREAM_Scale_verify(h_b, h_c, scalar, N);
644403adfb6SMatthew G Knepley   if (errorSTREAMkernel) {
645403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Scale:\t\tError detected in device STREAM_Scale, exiting\n");CHKERRQ(ierr);
646403adfb6SMatthew G Knepley     exit(-4000);
647403adfb6SMatthew G Knepley   } else {
648403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Scale:\t\tPass\n");CHKERRQ(ierr);
649403adfb6SMatthew G Knepley   }
650403adfb6SMatthew G Knepley 
651403adfb6SMatthew G Knepley   /* Initialize memory on the device */
652403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_a, 2.f, N);
653403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N);
654403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N);
655403adfb6SMatthew G Knepley 
656403adfb6SMatthew G Knepley   STREAM_Add<<<dimGrid,dimBlock>>>(d_a, d_b, d_c, N);
657403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_a, d_a, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
658403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_b, d_b, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
659403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_c, d_c, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
660403adfb6SMatthew G Knepley   errorSTREAMkernel = STREAM_Add_verify(h_a, h_b, h_c, N);
661403adfb6SMatthew G Knepley   if (errorSTREAMkernel) {
662403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Add:\t\tError detected in device STREAM_Add, exiting\n");CHKERRQ(ierr);
663403adfb6SMatthew G Knepley     exit(-5000);
664403adfb6SMatthew G Knepley   } else {
665403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Add:\t\tPass\n");CHKERRQ(ierr);
666403adfb6SMatthew G Knepley   }
667403adfb6SMatthew G Knepley 
668403adfb6SMatthew G Knepley   /* Initialize memory on the device */
669403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_a, 2.f, N);
670403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_b, .5f, N);
671403adfb6SMatthew G Knepley   set_array<<<dimGrid,dimBlock>>>(d_c, .5f, N);
672403adfb6SMatthew G Knepley 
673403adfb6SMatthew G Knepley   STREAM_Triad<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar, N);
674403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_a, d_a, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
675403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_b, d_b, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
676403adfb6SMatthew G Knepley   ierr              = cudaMemcpy(h_c, d_c, sizeof(float) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
677403adfb6SMatthew G Knepley   errorSTREAMkernel = STREAM_Triad_verify(h_b, h_c, h_a, scalar, N);
678403adfb6SMatthew G Knepley   if (errorSTREAMkernel) {
679403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Triad:\t\tError detected in device STREAM_Triad, exiting\n");CHKERRQ(ierr);
680403adfb6SMatthew G Knepley     exit(-6000);
681403adfb6SMatthew G Knepley   } else {
682403adfb6SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Triad:\t\tPass\n");CHKERRQ(ierr);
683403adfb6SMatthew G Knepley   }
684403adfb6SMatthew G Knepley 
685403adfb6SMatthew G Knepley   /* continue from here */
686403adfb6SMatthew G Knepley   printResultsReadable(times);
687403adfb6SMatthew G Knepley 
688403adfb6SMatthew G Knepley   //clean up timers
689403adfb6SMatthew G Knepley   ierr = cudaEventDestroy(stop);CHKERRQ(ierr);
690403adfb6SMatthew G Knepley   ierr = cudaEventDestroy(start);CHKERRQ(ierr);
691403adfb6SMatthew G Knepley 
692403adfb6SMatthew G Knepley   /* Free memory on device */
693403adfb6SMatthew G Knepley   ierr = cudaFree(d_a);CHKERRQ(ierr);
694403adfb6SMatthew G Knepley   ierr = cudaFree(d_b);CHKERRQ(ierr);
695403adfb6SMatthew G Knepley   ierr = cudaFree(d_c);CHKERRQ(ierr);
696403adfb6SMatthew G Knepley   PetscFunctionReturn(0);
697403adfb6SMatthew G Knepley }
698403adfb6SMatthew G Knepley 
699caccb7e3SMatthew G Knepley PetscErrorCode runStreamDouble(const PetscInt iNumThreadsPerBlock, PetscBool bDontUseGPUTiming)
700caccb7e3SMatthew G Knepley {
701caccb7e3SMatthew G Knepley   double         *d_a, *d_b, *d_c;
702caccb7e3SMatthew G Knepley   int            k;
703caccb7e3SMatthew G Knepley   float          times[8][NTIMES];
704caccb7e3SMatthew G Knepley   double         scalar;
705caccb7e3SMatthew G Knepley   PetscErrorCode ierr;
706caccb7e3SMatthew G Knepley 
707caccb7e3SMatthew G Knepley   PetscFunctionBegin;
708caccb7e3SMatthew G Knepley   /* Allocate memory on device */
709caccb7e3SMatthew G Knepley   ierr = cudaMalloc((void**)&d_a, sizeof(double)*N);CHKERRQ(ierr);
710caccb7e3SMatthew G Knepley   ierr = cudaMalloc((void**)&d_b, sizeof(double)*N);CHKERRQ(ierr);
711caccb7e3SMatthew G Knepley   ierr = cudaMalloc((void**)&d_c, sizeof(double)*N);CHKERRQ(ierr);
712caccb7e3SMatthew G Knepley 
713caccb7e3SMatthew G Knepley   /* Compute execution configuration */
714caccb7e3SMatthew G Knepley 
715caccb7e3SMatthew G Knepley   dim3 dimBlock(iNumThreadsPerBlock); /* (iNumThreadsPerBlock,1,1) */
716caccb7e3SMatthew G Knepley   dim3 dimGrid(N/dimBlock.x); /* (N/dimBlock.x,1,1) */
717caccb7e3SMatthew G Knepley   if (N % dimBlock.x != 0) dimGrid.x+=1;
718caccb7e3SMatthew G Knepley 
719caccb7e3SMatthew G Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, " Array size (double precision) = %u\n",N);CHKERRQ(ierr);
720caccb7e3SMatthew G Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, " using %u threads per block, %u blocks\n",dimBlock.x,dimGrid.x);CHKERRQ(ierr);
721caccb7e3SMatthew G Knepley 
722caccb7e3SMatthew G Knepley   /* Initialize memory on the device */
723caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_a, 2., N);
724caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N);
725caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N);
726caccb7e3SMatthew G Knepley 
727caccb7e3SMatthew G Knepley   /* --- MAIN LOOP --- repeat test cases NTIMES times --- */
728caccb7e3SMatthew G Knepley   PetscLogDouble cpuTimer = 0.0;
729caccb7e3SMatthew G Knepley   cudaEvent_t    start, stop;
730caccb7e3SMatthew G Knepley 
731caccb7e3SMatthew G Knepley   /* both timers report msec */
732caccb7e3SMatthew G Knepley   ierr = cudaEventCreate(&start);CHKERRQ(ierr); /* gpu timer facility */
733caccb7e3SMatthew G Knepley   ierr = cudaEventCreate(&stop);CHKERRQ(ierr);  /* gpu timer facility */
734caccb7e3SMatthew G Knepley 
735caccb7e3SMatthew G Knepley   scalar=3.0;
736caccb7e3SMatthew G Knepley   for (k = 0; k < NTIMES; ++k) {
737*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
738caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
739caccb7e3SMatthew G Knepley     STREAM_Copy_double<<<dimGrid,dimBlock>>>(d_a, d_c, N);
740caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
741caccb7e3SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
742caccb7e3SMatthew G Knepley     //get the the total elapsed time in ms
743caccb7e3SMatthew G Knepley     if (bDontUseGPUTiming) {
744*8563dfccSBarry Smith       PetscTimeAdd(&cpuTimer);
745caccb7e3SMatthew G Knepley       times[0][k] = cpuTimer;
746caccb7e3SMatthew G Knepley     } else {
747caccb7e3SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[0][k], start, stop);CHKERRQ(ierr);
748caccb7e3SMatthew G Knepley     }
749caccb7e3SMatthew G Knepley 
750caccb7e3SMatthew G Knepley     cpuTimer = 0.0;
751*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
752caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
753caccb7e3SMatthew G Knepley     STREAM_Copy_Optimized_double<<<dimGrid,dimBlock>>>(d_a, d_c, N);
754caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
755caccb7e3SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
756caccb7e3SMatthew G Knepley     //get the the total elapsed time in ms
757caccb7e3SMatthew G Knepley     if (bDontUseGPUTiming) {
758*8563dfccSBarry Smith       PetscTimeAdd(&cpuTimer);
759caccb7e3SMatthew G Knepley       times[1][k] = cpuTimer;
760caccb7e3SMatthew G Knepley     } else {
761caccb7e3SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[1][k], start, stop);CHKERRQ(ierr);
762caccb7e3SMatthew G Knepley     }
763caccb7e3SMatthew G Knepley 
764caccb7e3SMatthew G Knepley     cpuTimer = 0.0;
765*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
766caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
767caccb7e3SMatthew G Knepley     STREAM_Scale_double<<<dimGrid,dimBlock>>>(d_b, d_c, scalar,  N);
768caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
769caccb7e3SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
770caccb7e3SMatthew G Knepley     //get the the total elapsed time in ms
771*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
7726f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[2][k] = cpuTimer;
7736f2b61bcSKarl Rupp     else {
774caccb7e3SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[2][k], start, stop);CHKERRQ(ierr);
775caccb7e3SMatthew G Knepley     }
776caccb7e3SMatthew G Knepley 
777caccb7e3SMatthew G Knepley     cpuTimer = 0.0;
778*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
779caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
780caccb7e3SMatthew G Knepley     STREAM_Scale_Optimized_double<<<dimGrid,dimBlock>>>(d_b, d_c, scalar,  N);
781caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
782caccb7e3SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
783caccb7e3SMatthew G Knepley     //get the the total elapsed time in ms
784*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
7856f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[3][k] = cpuTimer;
7866f2b61bcSKarl Rupp     else {
787caccb7e3SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[2][k], start, stop);CHKERRQ(ierr);
788caccb7e3SMatthew G Knepley     }
789caccb7e3SMatthew G Knepley 
790caccb7e3SMatthew G Knepley     cpuTimer = 0.0;
791*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
792caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
793caccb7e3SMatthew G Knepley     STREAM_Add_double<<<dimGrid,dimBlock>>>(d_a, d_b, d_c,  N);
794caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
795caccb7e3SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
796caccb7e3SMatthew G Knepley     //get the the total elapsed time in ms
797*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
7986f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[4][k] = cpuTimer;
7996f2b61bcSKarl Rupp     else {
800caccb7e3SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[3][k], start, stop);CHKERRQ(ierr);
801caccb7e3SMatthew G Knepley     }
802caccb7e3SMatthew G Knepley 
803caccb7e3SMatthew G Knepley     cpuTimer = 0.0;
804*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
805caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
806caccb7e3SMatthew G Knepley     STREAM_Add_Optimized_double<<<dimGrid,dimBlock>>>(d_a, d_b, d_c,  N);
807caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
808caccb7e3SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
809caccb7e3SMatthew G Knepley     //get the the total elapsed time in ms
810*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
8116f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[5][k] = cpuTimer;
8126f2b61bcSKarl Rupp     else {
813caccb7e3SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[3][k], start, stop);CHKERRQ(ierr);
814caccb7e3SMatthew G Knepley     }
815caccb7e3SMatthew G Knepley 
816caccb7e3SMatthew G Knepley     cpuTimer = 0.0;
817*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
818caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
819caccb7e3SMatthew G Knepley     STREAM_Triad_double<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar,  N);
820caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
821caccb7e3SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
822caccb7e3SMatthew G Knepley     //get the the total elapsed time in ms
823*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
8246f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[6][k] = cpuTimer;
8256f2b61bcSKarl Rupp     else {
826caccb7e3SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[4][k], start, stop);CHKERRQ(ierr);
827caccb7e3SMatthew G Knepley     }
828caccb7e3SMatthew G Knepley 
829caccb7e3SMatthew G Knepley     cpuTimer = 0.0;
830*8563dfccSBarry Smith     PetscTimeSubtract(&cpuTimer);
831caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(start, 0);CHKERRQ(ierr);
832caccb7e3SMatthew G Knepley     STREAM_Triad_Optimized_double<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar,  N);
833caccb7e3SMatthew G Knepley     ierr = cudaEventRecord(stop, 0);CHKERRQ(ierr);
834caccb7e3SMatthew G Knepley     ierr = cudaEventSynchronize(stop);CHKERRQ(ierr);
835caccb7e3SMatthew G Knepley     //get the the total elapsed time in ms
836*8563dfccSBarry Smith     PetscTimeAdd(&cpuTimer);
8376f2b61bcSKarl Rupp     if (bDontUseGPUTiming) times[7][k] = cpuTimer;
8386f2b61bcSKarl Rupp     else {
839caccb7e3SMatthew G Knepley       ierr = cudaEventElapsedTime(&times[4][k], start, stop);CHKERRQ(ierr);
840caccb7e3SMatthew G Knepley     }
841caccb7e3SMatthew G Knepley 
842caccb7e3SMatthew G Knepley   }
843caccb7e3SMatthew G Knepley 
844caccb7e3SMatthew G Knepley   /* verify kernels */
845caccb7e3SMatthew G Knepley   double *h_a, *h_b, *h_c;
846caccb7e3SMatthew G Knepley   bool   errorSTREAMkernel = true;
847caccb7e3SMatthew G Knepley 
848caccb7e3SMatthew G Knepley   if ((h_a = (double*)calloc(N, sizeof(double))) == (double*)NULL) {
849caccb7e3SMatthew G Knepley     printf("Unable to allocate array h_a, exiting ...\n");
850caccb7e3SMatthew G Knepley     exit(1);
851caccb7e3SMatthew G Knepley   }
852caccb7e3SMatthew G Knepley   if ((h_b = (double*)calloc(N, sizeof(double))) == (double*)NULL) {
853caccb7e3SMatthew G Knepley     printf("Unable to allocate array h_b, exiting ...\n");
854caccb7e3SMatthew G Knepley     exit(1);
855caccb7e3SMatthew G Knepley   }
856caccb7e3SMatthew G Knepley 
857caccb7e3SMatthew G Knepley   if ((h_c = (double*)calloc(N, sizeof(double))) == (double*)NULL) {
858caccb7e3SMatthew G Knepley     printf("Unalbe to allocate array h_c, exiting ...\n");
859caccb7e3SMatthew G Knepley     exit(1);
860caccb7e3SMatthew G Knepley   }
861caccb7e3SMatthew G Knepley 
862caccb7e3SMatthew G Knepley   /*
863caccb7e3SMatthew G Knepley    * perform kernel, copy device memory into host memory and verify each
864caccb7e3SMatthew G Knepley    * device kernel output
865caccb7e3SMatthew G Knepley    */
866caccb7e3SMatthew G Knepley 
867caccb7e3SMatthew G Knepley   /* Initialize memory on the device */
868caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_a, 2., N);
869caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N);
870caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N);
871caccb7e3SMatthew G Knepley 
872caccb7e3SMatthew G Knepley   STREAM_Copy_double<<<dimGrid,dimBlock>>>(d_a, d_c, N);
873caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_a, d_a, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
874caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_c, d_c, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
875caccb7e3SMatthew G Knepley   errorSTREAMkernel = STREAM_Copy_verify_double(h_a, h_c, N);
876caccb7e3SMatthew G Knepley   if (errorSTREAMkernel) {
877caccb7e3SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy:\t\tError detected in device STREAM_Copy, exiting\n");CHKERRQ(ierr);
878caccb7e3SMatthew G Knepley     exit(-2000);
879caccb7e3SMatthew G Knepley   } else {
880caccb7e3SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy:\t\tPass\n");CHKERRQ(ierr);
881caccb7e3SMatthew G Knepley   }
882caccb7e3SMatthew G Knepley 
883caccb7e3SMatthew G Knepley   /* Initialize memory on the device */
884caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_a, 2., N);
885caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N);
886caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N);
887caccb7e3SMatthew G Knepley 
888caccb7e3SMatthew G Knepley   STREAM_Copy_Optimized_double<<<dimGrid,dimBlock>>>(d_a, d_c, N);
889caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_a, d_a, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
890caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_c, d_c, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
891caccb7e3SMatthew G Knepley   errorSTREAMkernel = STREAM_Copy_verify_double(h_a, h_c, N);
892caccb7e3SMatthew G Knepley   if (errorSTREAMkernel) {
893caccb7e3SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy_Optimized:\tError detected in device STREAM_Copy_Optimized, exiting\n");CHKERRQ(ierr);
894caccb7e3SMatthew G Knepley     exit(-3000);
895caccb7e3SMatthew G Knepley   } else {
896caccb7e3SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Copy_Optimized:\tPass\n");CHKERRQ(ierr);
897caccb7e3SMatthew G Knepley   }
898caccb7e3SMatthew G Knepley 
899caccb7e3SMatthew G Knepley   /* Initialize memory on the device */
900caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N);
901caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N);
902caccb7e3SMatthew G Knepley 
903caccb7e3SMatthew G Knepley   STREAM_Scale_double<<<dimGrid,dimBlock>>>(d_b, d_c, scalar, N);
904caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_b, d_b, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
905caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_c, d_c, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
906caccb7e3SMatthew G Knepley   errorSTREAMkernel = STREAM_Scale_verify_double(h_b, h_c, scalar, N);
907caccb7e3SMatthew G Knepley   if (errorSTREAMkernel) {
908caccb7e3SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Scale:\t\tError detected in device STREAM_Scale, exiting\n");CHKERRQ(ierr);
909caccb7e3SMatthew G Knepley     exit(-4000);
910caccb7e3SMatthew G Knepley   } else {
911caccb7e3SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Scale:\t\tPass\n");CHKERRQ(ierr);
912caccb7e3SMatthew G Knepley   }
913caccb7e3SMatthew G Knepley 
914caccb7e3SMatthew G Knepley   /* Initialize memory on the device */
915caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_a, 2., N);
916caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N);
917caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N);
918caccb7e3SMatthew G Knepley 
919caccb7e3SMatthew G Knepley   STREAM_Add_double<<<dimGrid,dimBlock>>>(d_a, d_b, d_c, N);
920caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_a, d_a, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
921caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_b, d_b, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
922caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_c, d_c, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
923caccb7e3SMatthew G Knepley   errorSTREAMkernel = STREAM_Add_verify_double(h_a, h_b, h_c, N);
924caccb7e3SMatthew G Knepley   if (errorSTREAMkernel) {
925caccb7e3SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Add:\t\tError detected in device STREAM_Add, exiting\n");CHKERRQ(ierr);
926caccb7e3SMatthew G Knepley     exit(-5000);
927caccb7e3SMatthew G Knepley   } else {
928caccb7e3SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Add:\t\tPass\n");CHKERRQ(ierr);
929caccb7e3SMatthew G Knepley   }
930caccb7e3SMatthew G Knepley 
931caccb7e3SMatthew G Knepley   /* Initialize memory on the device */
932caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_a, 2., N);
933caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_b, .5, N);
934caccb7e3SMatthew G Knepley   set_array_double<<<dimGrid,dimBlock>>>(d_c, .5, N);
935caccb7e3SMatthew G Knepley 
936caccb7e3SMatthew G Knepley   STREAM_Triad_double<<<dimGrid,dimBlock>>>(d_b, d_c, d_a, scalar, N);
937caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_a, d_a, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
938caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_b, d_b, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
939caccb7e3SMatthew G Knepley   ierr              = cudaMemcpy(h_c, d_c, sizeof(double) * N, cudaMemcpyDeviceToHost);CHKERRQ(ierr);
940caccb7e3SMatthew G Knepley   errorSTREAMkernel = STREAM_Triad_verify_double(h_b, h_c, h_a, scalar, N);
941caccb7e3SMatthew G Knepley   if (errorSTREAMkernel) {
942caccb7e3SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Triad:\t\tError detected in device STREAM_Triad, exiting\n");CHKERRQ(ierr);
943caccb7e3SMatthew G Knepley     exit(-6000);
944caccb7e3SMatthew G Knepley   } else {
945caccb7e3SMatthew G Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, " device STREAM_Triad:\t\tPass\n");CHKERRQ(ierr);
946caccb7e3SMatthew G Knepley   }
947caccb7e3SMatthew G Knepley 
948caccb7e3SMatthew G Knepley   /* continue from here */
949caccb7e3SMatthew G Knepley   printResultsReadable(times);
950caccb7e3SMatthew G Knepley 
951caccb7e3SMatthew G Knepley   //clean up timers
952caccb7e3SMatthew G Knepley   ierr = cudaEventDestroy(stop);CHKERRQ(ierr);
953caccb7e3SMatthew G Knepley   ierr = cudaEventDestroy(start);CHKERRQ(ierr);
954caccb7e3SMatthew G Knepley 
955caccb7e3SMatthew G Knepley   /* Free memory on device */
956caccb7e3SMatthew G Knepley   ierr = cudaFree(d_a);CHKERRQ(ierr);
957caccb7e3SMatthew G Knepley   ierr = cudaFree(d_b);CHKERRQ(ierr);
958caccb7e3SMatthew G Knepley   ierr = cudaFree(d_c);CHKERRQ(ierr);
959caccb7e3SMatthew G Knepley   PetscFunctionReturn(0);
960caccb7e3SMatthew G Knepley }
961caccb7e3SMatthew G Knepley 
962403adfb6SMatthew G Knepley ///////////////////////////////////////////////////////////////////////////
963403adfb6SMatthew G Knepley //Print Results to Screen and File
964403adfb6SMatthew G Knepley ///////////////////////////////////////////////////////////////////////////
965a6dfd86eSKarl Rupp PetscErrorCode printResultsReadable(float times[][NTIMES])
966a6dfd86eSKarl Rupp {
967403adfb6SMatthew G Knepley   PetscErrorCode ierr;
968403adfb6SMatthew G Knepley   PetscInt       j, k;
969caccb7e3SMatthew G Knepley   float          avgtime[8]          = {0., 0., 0., 0., 0., 0., 0., 0.};
970caccb7e3SMatthew G Knepley   float          maxtime[8]          = {0., 0., 0., 0., 0., 0., 0., 0.};
971caccb7e3SMatthew G Knepley   float          mintime[8]          = {1e30,1e30,1e30,1e30,1e30,1e30,1e30,1e30};
972caccb7e3SMatthew G Knepley   char           *label[8]           = {"Copy:      ", "Copy Opt.: ", "Scale:     ", "Scale Opt: ", "Add:       ", "Add Opt:   ", "Triad:     ", "Triad Opt: "};
973caccb7e3SMatthew G Knepley   float          bytes_per_kernel[8] = {
974403adfb6SMatthew G Knepley     2. * sizeof(float) * N,
975403adfb6SMatthew G Knepley     2. * sizeof(float) * N,
976403adfb6SMatthew G Knepley     2. * sizeof(float) * N,
977caccb7e3SMatthew G Knepley     2. * sizeof(float) * N,
978caccb7e3SMatthew G Knepley     3. * sizeof(float) * N,
979caccb7e3SMatthew G Knepley     3. * sizeof(float) * N,
980403adfb6SMatthew G Knepley     3. * sizeof(float) * N,
981403adfb6SMatthew G Knepley     3. * sizeof(float) * N
982403adfb6SMatthew G Knepley   };
983403adfb6SMatthew G Knepley 
984403adfb6SMatthew G Knepley   PetscFunctionBegin;
985403adfb6SMatthew G Knepley   /* --- SUMMARY --- */
9866f2b61bcSKarl Rupp   for (k = 1; k < NTIMES; ++k)   /* note -- skip first iteration */
987caccb7e3SMatthew G Knepley     for (j = 0; j < 8; ++j) {
988403adfb6SMatthew G Knepley       avgtime[j] = avgtime[j] + (1.e-03f * times[j][k]);
989403adfb6SMatthew G Knepley       mintime[j] = MIN(mintime[j], (1.e-03f * times[j][k]));
990403adfb6SMatthew G Knepley       maxtime[j] = MAX(maxtime[j], (1.e-03f * times[j][k]));
991403adfb6SMatthew G Knepley     }
992403adfb6SMatthew G Knepley 
993403adfb6SMatthew G Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, "Function    Rate (MB/s)    Avg time      Min time      Max time\n");CHKERRQ(ierr);
994403adfb6SMatthew G Knepley 
995caccb7e3SMatthew G Knepley   for (j = 0; j < 8; ++j) {
996403adfb6SMatthew G Knepley     avgtime[j] = avgtime[j]/(float)(NTIMES-1);
997403adfb6SMatthew G Knepley     ierr       = PetscPrintf(PETSC_COMM_SELF, "%s%11.4f  %11.6f  %12.6f  %12.6f\n", label[j], 1.0E-06 * bytes_per_kernel[j]/mintime[j], avgtime[j], mintime[j], maxtime[j]);CHKERRQ(ierr);
998403adfb6SMatthew G Knepley   }
999403adfb6SMatthew G Knepley   PetscFunctionReturn(0);
1000403adfb6SMatthew G Knepley }
1001