Actual source code: ex2.c

petsc-3.4.5 2014-06-29
  1: static const char help[] = "STREAM benchmark for pthread implemenentations\n\n";

  3: /*-----------------------------------------------------------------------*/
  4: /* Program: Stream                                                       */
  5: /* Revision: $Id: stream.c,v 5.9 2009/04/11 16:35:00 mccalpin Exp mccalpin $ */
  6: /* Original code developed by John D. McCalpin                           */
  7: /* Programmers: John D. McCalpin                                         */
  8: /*              Joe R. Zagar                                             */
  9: /*                                                                       */
 10: /* This program measures memory transfer rates in MB/s for simple        */
 11: /* computational kernels coded in C.                                     */
 12: /*-----------------------------------------------------------------------*/
 13: /* Copyright 1991-2005: John D. McCalpin                                 */
 14: /*-----------------------------------------------------------------------*/
 15: /* License:                                                              */
 16: /*  1. You are free to use this program and/or to redistribute           */
 17: /*     this program.                                                     */
 18: /*  2. You are free to modify this program for your own use,             */
 19: /*     including commercial use, subject to the publication              */
 20: /*     restrictions in item 3.                                           */
 21: /*  3. You are free to publish results obtained from running this        */
 22: /*     program, or from works that you derive from this program,         */
 23: /*     with the following limitations:                                   */
 24: /*     3a. In order to be referred to as "STREAM benchmark results",     */
 25: /*         published results must be in conformance to the STREAM        */
 26: /*         Run Rules, (briefly reviewed below) published at              */
 27: /*         http://www.cs.virginia.edu/stream/ref.html                    */
 28: /*         and incorporated herein by reference.                         */
 29: /*         As the copyright holder, John McCalpin retains the            */
 30: /*         right to determine conformity with the Run Rules.             */
 31: /*     3b. Results based on modified source code or on runs not in       */
 32: /*         accordance with the STREAM Run Rules must be clearly          */
 33: /*         labelled whenever they are published.  Examples of            */
 34: /*         proper labelling include:                                     */
 35: /*         "tuned STREAM benchmark results"                              */
 36: /*         "based on a variant of the STREAM benchmark code"             */
 37: /*         Other comparable, clear and reasonable labelling is           */
 38: /*         acceptable.                                                   */
 39: /*     3c. Submission of results to the STREAM benchmark web site        */
 40: /*         is encouraged, but not required.                              */
 41: /*  4. Use of this program or creation of derived works based on this    */
 42: /*     program constitutes acceptance of these licensing restrictions.   */
 43: /*  5. Absolutely no warranty is expressed or implied.                   */
 44: /*-----------------------------------------------------------------------*/
 45: # include <float.h>
 46: # include <sys/time.h>
 47: # include <petscsys.h>
 48: # include <petscthreadcomm.h>
 49: /* INSTRUCTIONS:
 50:  *
 51:  *      1) Stream requires a good bit of memory to run.  Adjust the
 52:  *          value of 'N' (below) to give a 'timing calibration' of
 53:  *          at least 20 clock-ticks.  This will provide rate estimates
 54:  *          that should be good to about 5% precision.
 55:  */

 57: #if !defined(N)
 58: #   define N  2000000
 59: #endif
 60: #if !defined(NTIMES)
 61: #   define NTIMES       50
 62: #endif
 63: #if !defined(OFFSET)
 64: #   define OFFSET       0
 65: #endif

 67: /*
 68:  *      3) Compile the code with full optimization.  Many compilers
 69:  *         generate unreasonably bad code before the optimizer tightens
 70:  *         things up.  If the results are unreasonably good, on the
 71:  *         other hand, the optimizer might be too smart for me!
 72:  *
 73:  *         Try compiling with:
 74:  *               cc -O stream_omp.c -o stream_omp
 75:  *
 76:  *         This is known to work on Cray, SGI, IBM, and Sun machines.
 77:  *
 78:  *
 79:  *      4) Mail the results to mccalpin@cs.virginia.edu
 80:  *         Be sure to include:
 81:  *              a) computer hardware model number and software revision
 82:  *              b) the compiler flags
 83:  *              c) all of the output from the test case.
 84:  * Thanks!
 85:  *
 86:  */

 88: # define HLINE "-------------------------------------------------------------\n"

 90: # if !defined(MIN)
 91: # define MIN(x,y) ((x)<(y) ? (x) : (y))
 92: # endif
 93: # if !defined(MAX)
 94: # define MAX(x,y) ((x)>(y) ? (x) : (y))
 95: # endif

 97: #if !defined(STATIC_ALLOC)
 98: #  define STATIC_ALLOC 1
 99: #endif

101: #if STATIC_ALLOC
102: static double a[N+OFFSET],
103:               b[N+OFFSET],
104:               c[N+OFFSET];
105: #endif

107: static double avgtime[4] = {0}, maxtime[4] = {0},
108:               mintime[4] = {FLT_MAX,FLT_MAX,FLT_MAX,FLT_MAX};

110: static const char *label[4] = {"Copy:      ", "Scale:     ","Add:       ", "Triad:     "};

112: static double bytes[4] = {
113:   2 * sizeof(double) * N,
114:   2 * sizeof(double) * N,
115:   3 * sizeof(double) * N,
116:   3 * sizeof(double) * N
117: };

119: double mysecond();
120: void checkSTREAMresults();

122: #if !defined(WITH_PTHREADS)
123: #  define WITH_PTHREADS 1
124: #endif

126: #if !STATIC_ALLOC
127: double *a, *b, *c;
128: #endif

130: #if WITH_PTHREADS
131: PetscInt nworkThreads;
132: PetscInt *trstarts;

134: PetscErrorCode tuned_STREAM_2A_Kernel(PetscInt myrank)
135: {
136:   int i;

138:   for (i=trstarts[myrank]; i<trstarts[myrank+1]; i++) a[i] = 2.0E0*a[i];

140:   return(0);
141: }

143: PetscErrorCode tuned_STREAM_Initialize_Kernel(PetscInt myrank)
144: {
145:   int i;

147:   for (i=trstarts[myrank]; i<trstarts[myrank+1]; i++) {
148:     a[i] = 1.0;
149:     b[i] = 2.0;
150:     c[i] = 0.0;
151:   }
152:   return(0);
153: }

155: PetscErrorCode tuned_STREAM_Copy_Kernel(PetscInt myrank)
156: {
157:   int j;

159:   for (j=trstarts[myrank]; j<trstarts[myrank+1]; j++) c[j] = a[j];
160:   return(0);
161: }

163: PetscErrorCode tuned_STREAM_Scale_Kernel(PetscInt myrank,double *scalarp)
164: {
165:   double scalar = *scalarp;
166:   int    j;

168:   for (j=trstarts[myrank]; j<trstarts[myrank+1]; j++) b[j] = scalar*c[j];
169:   return(0);
170: }

172: PetscErrorCode tuned_STREAM_Add_Kernel(PetscInt myrank)
173: {
174:   int j;

176:   for (j=trstarts[myrank]; j<trstarts[myrank+1]; j++) c[j] = a[j]+b[j];

178:   return(0);
179: }

181: PetscErrorCode tuned_STREAM_Triad_Kernel(PetscInt myrank,double *scalarp)
182: {
183:   double scalar = *scalarp;
184:   int    j;

186:   for (j=trstarts[myrank]; j<trstarts[myrank+1]; j++) a[j] = b[j]+scalar*c[j];

188:   return(0);
189: }
190: #endif

192: int main(int argc,char *argv[])
193: {
195:   int            quantum, checktick();
196:   int            BytesPerWord;
197:   int            j, k;
198:   double         scalar=3.0, t, times[4][NTIMES];

200:   PetscInitialize(&argc,&argv,0,help);
201:   /* --- SETUP --- determine precision and check timing --- */

203:   /*printf(HLINE);
204:     printf("STREAM version $Revision: 5.9 $\n");
205:     printf(HLINE); */
206:   BytesPerWord = sizeof(double);
207:   printf("This system uses %d bytes per DOUBLE PRECISION word.\n",BytesPerWord);

209:   printf(HLINE);
210: #if defined(NO_LONG_LONG)
211:   printf("Array size = %d, Offset = %d\n", N, OFFSET);
212: #else
213:   printf("Array size = %llu, Offset = %d\n", (unsigned long long) N, OFFSET);
214: #endif

216:   printf("Total memory required = %.1f MB.\n",(3.0 * BytesPerWord) * ((double)N / 1048576.0));
217:   printf("Each test is run %d times, but only\n", NTIMES);
218:   printf("the *best* time for each is used.\n");

220:   printf(HLINE);

222: #if !STATIC_ALLOC
223:   a = malloc((N+OFFSET)*sizeof(double));
224:   b = malloc((N+OFFSET)*sizeof(double));
225:   c = malloc((N+OFFSET)*sizeof(double));
226: #endif

228: #if WITH_PTHREADS
229:   PetscThreadCommGetNThreads(PETSC_COMM_WORLD,&nworkThreads);
230:   PetscMalloc((nworkThreads+1)*sizeof(PetscInt),&trstarts);
231:   PetscInt  Q,R,nloc;
232:   PetscBool S;
233:   Q           = (N+OFFSET)/nworkThreads;
234:   R           = (N+OFFSET) - Q*nworkThreads;
235:   trstarts[0] = 0;
236:   for (j=0; j < nworkThreads; j++) {
237:     S             = (PetscBool)(j < R);
238:     nloc          = S ? Q+1 : Q;
239:     trstarts[j+1] = trstarts[j]+nloc;
240:   }
241:   PetscThreadCommRunKernel(PETSC_COMM_WORLD,(PetscThreadKernel)tuned_STREAM_Initialize_Kernel,1,&scalar);
242:   PetscThreadCommBarrier(PETSC_COMM_WORLD);
243: # else
244:   for (j=0; j<N; j++) {
245:     a[j] = 1.0;
246:     b[j] = 2.0;
247:     c[j] = 0.0;
248:   }
249: #endif

251:   /*printf(HLINE);*/

253:   /* Get initial value for system clock. */
254:   if  ((quantum = checktick()) >= 1) ;
255:   /*      printf("Your clock granularity/precision appears to be %d microseconds.\n", quantum); */
256:   else quantum = 1;
257:   /*   printf("Your clock granularity appears to be less than one microsecond.\n"); */

259:   t = mysecond();

261: #if WITH_PTHREADS
262:   PetscThreadCommRunKernel(PETSC_COMM_WORLD,(PetscThreadKernel)tuned_STREAM_2A_Kernel,0);
263:   PetscThreadCommBarrier(PETSC_COMM_WORLD);
264: #else
265:   for (j = 0; j < N; j++) a[j] = 2.0E0 * a[j];
266: #endif
267:   t = 1.0E6 * (mysecond() - t);

269:   /*    printf("Each test below will take on the order of %d microseconds.\n", (int)t);
270:   printf("   (= %d clock ticks)\n", (int) (t/quantum));
271:   printf("Increase the size of the arrays if this shows that\n");
272:   printf("you are not getting at least 20 clock ticks per test.\n");

274:   printf(HLINE);
275:   */
276:   /*  --- MAIN LOOP --- repeat test cases NTIMES times --- */

278:   for (k=0; k<NTIMES; k++) {
279:     times[0][k] = mysecond();
280: #if WITH_PTHREADS
281:     PetscThreadCommRunKernel(PETSC_COMM_WORLD,(PetscThreadKernel)tuned_STREAM_Copy_Kernel,0);
282:     PetscThreadCommBarrier(PETSC_COMM_WORLD);
283: #else
284:     for (j=0; j<N; j++) c[j] = a[j];
285: #endif
286:     times[0][k] = mysecond() - times[0][k];

288:     times[1][k] = mysecond();
289: #if WITH_PTHREADS
290:     PetscThreadCommRunKernel(PETSC_COMM_WORLD,(PetscThreadKernel)tuned_STREAM_Scale_Kernel,1,&scalar);
291:     PetscThreadCommBarrier(PETSC_COMM_WORLD);
292: #else
293:     for (j=0; j<N; j++) b[j] = scalar*c[j];
294: #endif
295:     times[1][k] = mysecond() - times[1][k];

297:     times[2][k] = mysecond();
298: #if WITH_PTHREADS
299:     PetscThreadCommRunKernel(PETSC_COMM_WORLD,(PetscThreadKernel)tuned_STREAM_Add_Kernel,0);
300:     PetscThreadCommBarrier(PETSC_COMM_WORLD);
301: #else
302:     for (j=0; j<N; j++) c[j] = a[j]+b[j];
303: #endif
304:     times[2][k] = mysecond() - times[2][k];

306:     times[3][k] = mysecond();
307: #if WITH_PTHREADS
308:     PetscThreadCommRunKernel(PETSC_COMM_WORLD,(PetscThreadKernel)tuned_STREAM_Triad_Kernel,1,&scalar);
309:     PetscThreadCommBarrier(PETSC_COMM_WORLD);
310: #else
311:     for (j=0; j<N; j++) a[j] = b[j]+scalar*c[j];
312: #endif
313:     times[3][k] = mysecond() - times[3][k];
314:   }

316:   /*  --- SUMMARY --- */

318:   for (k=1; k<NTIMES; k++)     /* note -- skip first iteration */
319:     for (j=0; j<4; j++) {
320:       avgtime[j] = avgtime[j] + times[j][k];
321:       mintime[j] = MIN(mintime[j], times[j][k]);
322:       maxtime[j] = MAX(maxtime[j], times[j][k]);
323:     }

325:   printf("Function      Rate (MB/s) \n");
326:   for (j=0; j<4; j++) {
327:     avgtime[j] = avgtime[j]/(double)(NTIMES-1);

329:     printf("%s%11.4f  \n", label[j], 1.0E-06 * bytes[j]/mintime[j]);
330:   }
331:   /* printf(HLINE);*/
332: #if WITH_PTHREADS
333:   PetscThreadCommBarrier(PETSC_COMM_WORLD);
334: #endif
335:   /* --- Check Results --- */
336:   checkSTREAMresults();
337:   /*    printf(HLINE);*/
338:   PetscFinalize();
339:   return 0;
340: }

342: # define        M        20

344: int checktick()
345: {
346:   int    i, minDelta, Delta;
347:   double t1, t2, timesfound[M];

349:   /*  Collect a sequence of M unique time values from the system. */

351:   for (i = 0; i < M; i++) {
352:     t1 = mysecond();
353:     while (((t2=mysecond()) - t1) < 1.0E-6) ;
354:     timesfound[i] = t1 = t2;
355:   }

357:   /*
358:    * Determine the minimum difference between these M values.
359:    * This result will be our estimate (in microseconds) for the
360:    * clock granularity.
361:    */

363:   minDelta = 1000000;
364:   for (i = 1; i < M; i++) {
365:     Delta    = (int)(1.0E6 * (timesfound[i]-timesfound[i-1]));
366:     minDelta = MIN(minDelta, MAX(Delta,0));
367:   }

369:   return(minDelta);
370: }

372: /* A gettimeofday routine to give access to the wall
373:    clock timer on most UNIX-like systems.  */

375: #include <sys/time.h>

377: double mysecond()
378: {
379:   struct timeval  tp;
380:   struct timezone tzp;
381:   int             i=0;

383:   i = gettimeofday(&tp,&tzp);
384:   return ((double) tp.tv_sec + (double) tp.tv_usec * 1.e-6);
385: }

387: void checkSTREAMresults()
388: {
389:   double aj,bj,cj,scalar;
390:   double asum,bsum,csum;
391:   double epsilon;
392:   int    j,k;

394:   /* reproduce initialization */
395:   aj = 1.0;
396:   bj = 2.0;
397:   cj = 0.0;
398:   /* a[] is modified during timing check */
399:   aj = 2.0E0 * aj;
400:   /* now execute timing loop */
401:   scalar = 3.0;
402:   for (k=0; k<NTIMES; k++) {
403:     cj = aj;
404:     bj = scalar*cj;
405:     cj = aj+bj;
406:     aj = bj+scalar*cj;
407:   }
408:   aj = aj * (double) (N);
409:   bj = bj * (double) (N);
410:   cj = cj * (double) (N);

412:   asum = 0.0;
413:   bsum = 0.0;
414:   csum = 0.0;
415:   for (j=0; j<N; j++) {
416:     asum += a[j];
417:     bsum += b[j];
418:     csum += c[j];
419:   }
420: #if defined(VERBOSE)
421:   printf ("Results Comparison: \n");
422:   printf ("        Expected  : %f %f %f \n",aj,bj,cj);
423:   printf ("        Observed  : %f %f %f \n",asum,bsum,csum);
424: #endif

426: #if !defined(abs)
427: #define abs(a) ((a) >= 0 ? (a) : -(a))
428: #endif
429:   epsilon = 1.e-8;

431:   if (abs(aj-asum)/asum > epsilon) {
432:     printf ("Failed Validation on array a[]\n");
433:     printf ("        Expected  : %f \n",aj);
434:     printf ("        Observed  : %f \n",asum);
435:   }
436:   if (abs(bj-bsum)/bsum > epsilon) {
437:     printf ("Failed Validation on array b[]\n");
438:     printf ("        Expected  : %f \n",bj);
439:     printf ("        Observed  : %f \n",bsum);
440:   }
441:   if (abs(cj-csum)/csum > epsilon) {
442:     printf ("Failed Validation on array c[]\n");
443:     printf ("        Expected  : %f \n",cj);
444:     printf ("        Observed  : %f \n",csum);
445:   } /* else printf ("Solution Validates\n"); */
446: }