Actual source code: ex2hip.hip.cpp

  1: static char help[] = "Benchmarking hipPointerGetAttributes() time\n";
  2: /*
  3:   Running example on Crusher at OLCF:
  4:     # run with 1 mpi rank (-n1), 32 CPUs (-c32), and map the process to CPU 0 and GPU 0
  5:   $ srun -n1 -c32 --cpu-bind=map_cpu:0 --gpus-per-node=8 --gpu-bind=map_gpu:0 ./ex2hip
  6:     Average hipPointerGetAttributes() time = 0.24 microseconds
  7: */
  8: #include <petscsys.h>
  9: #include <petscdevice_hip.h>

 11: int main(int argc, char **argv)
 12: {
 13:   PetscInt              i, n = 4000;
 14:   hipError_t            cerr;
 15:   PetscScalar         **ptrs;
 16:   PetscLogDouble        tstart, tend, time;
 17:   hipPointerAttribute_t attr;

 20:   PetscInitialize(&argc, &argv, (char *)0, help);
 21:   PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL);
 22:   hipStreamSynchronize(NULL); /* Initialize HIP runtime to get more accurate timing below */

 24:   PetscMalloc1(n, &ptrs);
 25:   for (i = 0; i < n; i++) {
 26:     if (i % 2) PetscMalloc1(i + 16, &ptrs[i]);
 27:     else hipMalloc((void **)&ptrs[i], (i + 16) * sizeof(PetscScalar));
 28:   }

 30:   PetscTime(&tstart);
 31:   for (i = 0; i < n; i++) {
 32:     cerr = hipPointerGetAttributes(&attr, ptrs[i]);
 33:     if (cerr) cerr = hipGetLastError();
 34:   }
 35:   PetscTime(&tend);
 36:   time = (tend - tstart) * 1e6 / n;

 38:   PetscPrintf(PETSC_COMM_WORLD, "Average hipPointerGetAttributes() time = %.2f microseconds\n", time);

 40:   for (i = 0; i < n; i++) {
 41:     if (i % 2) PetscFree(ptrs[i]);
 42:     else hipFree(ptrs[i]);
 43:   }
 44:   PetscFree(ptrs);
 45:   PetscFinalize();
 46:   return 0;
 47: }

 49: /*TEST
 50:   build:
 51:     requires: hip

 53:   test:
 54:     requires: hip
 55:     args: -n 2
 56:     output_file: output/empty.out
 57:     filter: grep "DOES_NOT_EXIST"

 59: TEST*/