Actual source code: performance.c
petsc-3.13.6 2020-09-29
1: static char help[] = "Time vector operations on GPU\n";
2: /* This program produces the results for Argonne Technical Report ANL-19/41.
3: The technical report and resources for generating data can be found in the
4: repository: https://gitlab.com/hannah_mairs/summit-performance */
6: #include <petscvec.h>
8: int main(int argc, char **argv)
9: {
11: Vec v,w,x;
12: PetscInt n=15;
13: PetscScalar val;
14: PetscReal norm1,norm2;
15: PetscRandom rctx;
16: PetscLogStage stage;
18: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
19: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
20: PetscRandomCreate(PETSC_COMM_WORLD,&rctx);
21: PetscRandomSetFromOptions(rctx);
22: VecCreate(PETSC_COMM_WORLD,&v);
23: VecSetSizes(v,PETSC_DECIDE,n);
24: VecSetFromOptions(v);
25: VecDuplicate(v,&w);
26: VecSetRandom(v,rctx);
27: VecSetRandom(w,rctx);
29: /* create dummy vector to clear cache */
30: VecCreate(PETSC_COMM_WORLD,&x);
31: VecSetSizes(x,PETSC_DECIDE,10000000);
32: VecSetFromOptions(x);
33: VecSetRandom(x,rctx);
35: /* send v to GPU */
36: PetscBarrier(NULL);
37: VecNorm(v,NORM_1,&norm1);
39: /* register a stage work on GPU */
40: PetscLogStageRegister("Work on GPU", &stage);
41: PetscLogStagePush(stage);
42: VecNorm(w,NORM_1,&norm1); /* send w to GPU */
43: VecNorm(x,NORM_1,&norm1); /* clear cache */
44: PetscBarrier(NULL);
45: VecAXPY(w,1.0,v);
46: VecNorm(x,NORM_INFINITY,&norm1);
47: PetscBarrier(NULL);
48: VecDot(w,v,&val);
49: VecNorm(x,NORM_1,&norm1);
50: PetscBarrier(NULL);
51: VecSet(v,0.0);
52: VecNorm(x,NORM_2,&norm2);
53: PetscBarrier(NULL);
54: VecCopy(v,w);
55: PetscLogStagePop();
57: PetscPrintf(PETSC_COMM_WORLD,"Test completed successfully!\n");
58: VecDestroy(&v);
59: VecDestroy(&w);
60: VecDestroy(&x);
61: PetscRandomDestroy(&rctx);
62: PetscFinalize();
63: return ierr;
64: }
66: /*TEST
68: test:
69: suffix: cuda
70: nsize: 2
71: args: -vec_type mpicuda
72: requires: cuda
74: TEST*/