1c4762a1bSJed Brown static char help[] = "Test scalability of PetscHSetIJ hash set.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscsys.h> 4c4762a1bSJed Brown #include <petsctime.h> 5c4762a1bSJed Brown #include <petsc/private/hashsetij.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown int main(int argc, char **argv) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown PetscHSetIJ table; 10c4762a1bSJed Brown PetscInt N = 0, i, j, n; 11c4762a1bSJed Brown PetscHashIJKey key; 12c4762a1bSJed Brown PetscBool flag; 13c4762a1bSJed Brown PetscLogDouble t_add = 0; 14c4762a1bSJed Brown PetscLogDouble t_has = 0; 15c4762a1bSJed Brown PetscLogDouble t_del = 0; 16c4762a1bSJed Brown PetscErrorCode ierr; 17c4762a1bSJed Brown 18c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIJCreate(&table)); 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* The following line silences warnings from Clang Static Analyzer */ 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIJResize(table,0)); 24c4762a1bSJed Brown 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTimeSubtract(&t_add)); 26c4762a1bSJed Brown for (i = 0; i < N; ++i) { 27c4762a1bSJed Brown for (j = 0; j < N; ++j) { 28c4762a1bSJed Brown key.i = PetscMin(i, j); 29c4762a1bSJed Brown key.j = PetscMax(i, j); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIJQueryAdd(table, key, &flag)); 31c4762a1bSJed Brown } 32c4762a1bSJed Brown } 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTimeAdd(&t_add)); 34c4762a1bSJed Brown 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIJGetSize(table,&n)); 36c4762a1bSJed Brown 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTimeSubtract(&t_has)); 38c4762a1bSJed Brown for (i = 0; i < N; ++i) { 39c4762a1bSJed Brown for (j = 0; j < N; ++j) { 40c4762a1bSJed Brown key.i = i; 41c4762a1bSJed Brown key.j = j; 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIJHas(table, key, &flag)); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown } 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTimeAdd(&t_has)); 46c4762a1bSJed Brown 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTimeSubtract(&t_del)); 48c4762a1bSJed Brown for (i = 0; i < N; ++i) { 49c4762a1bSJed Brown for (j = 0; j < N; ++j) { 50c4762a1bSJed Brown key.i = i; 51c4762a1bSJed Brown key.j = j; 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIJQueryDel(table, key, &flag)); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown } 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTimeAdd(&t_del)); 56c4762a1bSJed Brown 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"N = %" PetscInt_FMT " - table size: %" PetscInt_FMT ", add: %g, has: %g, del: %g\n",N,n,t_add,t_has,t_del)); 58c4762a1bSJed Brown 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIJDestroy(&table)); 60c4762a1bSJed Brown ierr = PetscFinalize(); 61c4762a1bSJed Brown return ierr; 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /*TEST 65c4762a1bSJed Brown 66c4762a1bSJed Brown test: 67c4762a1bSJed Brown args: -N 32 68c4762a1bSJed Brown 69c4762a1bSJed Brown TEST*/ 70