1c4762a1bSJed Brown static char help[] = "Test scalability of PetscHSetI hash set.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscsys.h> 4c4762a1bSJed Brown #include <petsctime.h> 5c4762a1bSJed Brown #include <petsc/private/hashseti.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown int main(int argc, char **argv) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown PetscHSetI table; 10c4762a1bSJed Brown PetscInt N = 0, i, j, n; 11c4762a1bSJed Brown PetscBool flag; 12c4762a1bSJed Brown PetscLogDouble t_add = 0; 13c4762a1bSJed Brown PetscLogDouble t_has = 0; 14c4762a1bSJed Brown PetscLogDouble t_del = 0; 15c4762a1bSJed Brown 16*327415f7SBarry Smith PetscFunctionBeginUser; 179566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 199566063dSJacob Faibussowitsch PetscCall(PetscHSetICreate(&table)); 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* The following line silences warnings from Clang Static Analyzer */ 229566063dSJacob Faibussowitsch PetscCall(PetscHSetIResize(table,0)); 23c4762a1bSJed Brown 249566063dSJacob Faibussowitsch PetscCall(PetscTimeSubtract(&t_add)); 25c4762a1bSJed Brown for (i = 0; i < N; ++i) { 26c4762a1bSJed Brown for (j = 0; j < N; ++j) { 27c4762a1bSJed Brown PetscInt key = i+j*N; 289566063dSJacob Faibussowitsch PetscCall(PetscHSetIQueryAdd(table, key, &flag)); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown } 319566063dSJacob Faibussowitsch PetscCall(PetscTimeAdd(&t_add)); 32c4762a1bSJed Brown 339566063dSJacob Faibussowitsch PetscCall(PetscHSetIGetSize(table,&n)); 34c4762a1bSJed Brown 359566063dSJacob Faibussowitsch PetscCall(PetscTimeSubtract(&t_has)); 36c4762a1bSJed Brown for (i = 0; i < N; ++i) { 37c4762a1bSJed Brown for (j = 0; j < N; ++j) { 38c4762a1bSJed Brown PetscInt key = i+j*N; 399566063dSJacob Faibussowitsch PetscCall(PetscHSetIHas(table, key, &flag)); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown } 429566063dSJacob Faibussowitsch PetscCall(PetscTimeAdd(&t_has)); 43c4762a1bSJed Brown 449566063dSJacob Faibussowitsch PetscCall(PetscTimeSubtract(&t_del)); 45c4762a1bSJed Brown for (i = 0; i < N; ++i) { 46c4762a1bSJed Brown for (j = 0; j < N; ++j) { 47c4762a1bSJed Brown PetscInt key = i+j*N; 489566063dSJacob Faibussowitsch PetscCall(PetscHSetIQueryDel(table, key, &flag)); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown } 519566063dSJacob Faibussowitsch PetscCall(PetscTimeAdd(&t_del)); 52c4762a1bSJed Brown 539566063dSJacob 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)); 54c4762a1bSJed Brown 559566063dSJacob Faibussowitsch PetscCall(PetscHSetIDestroy(&table)); 569566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 57b122ec5aSJacob Faibussowitsch return 0; 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown /*TEST 61c4762a1bSJed Brown 62c4762a1bSJed Brown test: 63c4762a1bSJed Brown args: -N 32 64c4762a1bSJed Brown 65c4762a1bSJed Brown TEST*/ 66