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 17*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 18*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 19*9566063dSJacob Faibussowitsch PetscCall(PetscHSetIJCreate(&table)); 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* The following line silences warnings from Clang Static Analyzer */ 22*9566063dSJacob Faibussowitsch PetscCall(PetscHSetIJResize(table,0)); 23c4762a1bSJed Brown 24*9566063dSJacob Faibussowitsch PetscCall(PetscTimeSubtract(&t_add)); 25c4762a1bSJed Brown for (i = 0; i < N; ++i) { 26c4762a1bSJed Brown for (j = 0; j < N; ++j) { 27c4762a1bSJed Brown key.i = PetscMin(i, j); 28c4762a1bSJed Brown key.j = PetscMax(i, j); 29*9566063dSJacob Faibussowitsch PetscCall(PetscHSetIJQueryAdd(table, key, &flag)); 30c4762a1bSJed Brown } 31c4762a1bSJed Brown } 32*9566063dSJacob Faibussowitsch PetscCall(PetscTimeAdd(&t_add)); 33c4762a1bSJed Brown 34*9566063dSJacob Faibussowitsch PetscCall(PetscHSetIJGetSize(table,&n)); 35c4762a1bSJed Brown 36*9566063dSJacob Faibussowitsch PetscCall(PetscTimeSubtract(&t_has)); 37c4762a1bSJed Brown for (i = 0; i < N; ++i) { 38c4762a1bSJed Brown for (j = 0; j < N; ++j) { 39c4762a1bSJed Brown key.i = i; 40c4762a1bSJed Brown key.j = j; 41*9566063dSJacob Faibussowitsch PetscCall(PetscHSetIJHas(table, key, &flag)); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown } 44*9566063dSJacob Faibussowitsch PetscCall(PetscTimeAdd(&t_has)); 45c4762a1bSJed Brown 46*9566063dSJacob Faibussowitsch PetscCall(PetscTimeSubtract(&t_del)); 47c4762a1bSJed Brown for (i = 0; i < N; ++i) { 48c4762a1bSJed Brown for (j = 0; j < N; ++j) { 49c4762a1bSJed Brown key.i = i; 50c4762a1bSJed Brown key.j = j; 51*9566063dSJacob Faibussowitsch PetscCall(PetscHSetIJQueryDel(table, key, &flag)); 52c4762a1bSJed Brown } 53c4762a1bSJed Brown } 54*9566063dSJacob Faibussowitsch PetscCall(PetscTimeAdd(&t_del)); 55c4762a1bSJed Brown 56*9566063dSJacob Faibussowitsch PetscCall(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)); 57c4762a1bSJed Brown 58*9566063dSJacob Faibussowitsch PetscCall(PetscHSetIJDestroy(&table)); 59*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 60b122ec5aSJacob Faibussowitsch return 0; 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 63c4762a1bSJed Brown /*TEST 64c4762a1bSJed Brown 65c4762a1bSJed Brown test: 66c4762a1bSJed Brown args: -N 32 67c4762a1bSJed Brown 68c4762a1bSJed Brown TEST*/ 69