xref: /petsc/src/benchmarks/PetscMalloc.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1d8e9fea7SSatish Balay 
2c6db04a5SJed Brown #include <petscsys.h>
38563dfccSBarry Smith #include <petsctime.h>
4173c0623SSatish Balay 
5173c0623SSatish Balay int main(int argc,char **argv)
6173c0623SSatish Balay {
7b0a32e0cSBarry Smith   PetscLogDouble x,y;
847794344SBarry Smith   double         value;
998d6af09SSatish Balay   void           *arr[1000],*dummy;
10a438ae71SBarry Smith   int            i,rand1[1000],rand2[1000];
11a438ae71SBarry Smith   PetscErrorCode ierr;
1277c4ece6SBarry Smith   PetscRandom    r;
13ace3abfcSBarry Smith   PetscBool      flg;
14173c0623SSatish Balay 
15a438ae71SBarry Smith   ierr = PetscInitialize(&argc,&argv,0,0);if (ierr) return ierr;
16*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&r));
17*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(r));
18173c0623SSatish Balay   for (i=0; i<1000; i++) {
19*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(r,&value));
20173c0623SSatish Balay     rand1[i] = (int)(value* 144327);
21*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(r,&value));
22173c0623SSatish Balay     rand2[i] = (int)(value* 144327);
23173c0623SSatish Balay   }
24173c0623SSatish Balay 
2598d6af09SSatish Balay   /* Take care of paging effects */
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(100,&dummy));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(dummy));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTime(&x));
2998d6af09SSatish Balay 
30173c0623SSatish Balay   /* Do all mallocs */
31173c0623SSatish Balay   for (i=0; i< 1000; i++) {
32*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(rand1[i],&arr[i]));
33173c0623SSatish Balay   }
34173c0623SSatish Balay 
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTime(&x));
36173c0623SSatish Balay 
37173c0623SSatish Balay   /* Do some frees */
38173c0623SSatish Balay   for (i=0; i< 1000; i+=2) {
39*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(arr[i]));
40173c0623SSatish Balay   }
41173c0623SSatish Balay 
42173c0623SSatish Balay   /* Do some mallocs */
43173c0623SSatish Balay   for (i=0; i< 1000; i+=2) {
44*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(rand2[i],&arr[i]));
45173c0623SSatish Balay   }
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTime(&y));
47173c0623SSatish Balay 
48173c0623SSatish Balay   for (i=0; i< 1000; i++) {
49*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(arr[i]));
50173c0623SSatish Balay   }
51173c0623SSatish Balay 
52b4d8b9abSSatish Balay   fprintf(stdout,"%-15s : %e sec, with options : ","PetscMalloc",(y-x)/500.0);
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,"-malloc",&flg));
548caf3d72SBarry Smith   if (flg) fprintf(stdout,"-malloc ");
55b4d8b9abSSatish Balay   fprintf(stdout,"\n");
56173c0623SSatish Balay 
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&r));
58f3fe499bSBarry Smith   ierr = PetscFinalize();
5926f47effSBarry Smith   return ierr;
60173c0623SSatish Balay }
61