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 PetscErrorCode ierr; 16c4762a1bSJed Brown 17c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 18*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetICreate(&table)); 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* The following line silences warnings from Clang Static Analyzer */ 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIResize(table,0)); 23c4762a1bSJed Brown 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIQueryAdd(table, key, &flag)); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown } 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTimeAdd(&t_add)); 32c4762a1bSJed Brown 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIGetSize(table,&n)); 34c4762a1bSJed Brown 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIHas(table, key, &flag)); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown } 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTimeAdd(&t_has)); 43c4762a1bSJed Brown 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIQueryDel(table, key, &flag)); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown } 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTimeAdd(&t_del)); 52c4762a1bSJed Brown 53*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)); 54c4762a1bSJed Brown 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHSetIDestroy(&table)); 56c4762a1bSJed Brown ierr = PetscFinalize(); 57c4762a1bSJed Brown return ierr; 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown /*TEST 61c4762a1bSJed Brown 62c4762a1bSJed Brown test: 63c4762a1bSJed Brown args: -N 32 64c4762a1bSJed Brown 65c4762a1bSJed Brown TEST*/ 66