xref: /petsc/src/benchmarks/PetscMalloc.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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];
1177c4ece6SBarry Smith   PetscRandom    r;
12ace3abfcSBarry Smith   PetscBool      flg;
13173c0623SSatish Balay 
14*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,0,0));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&r));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(r));
17173c0623SSatish Balay   for (i=0; i<1000; i++) {
185f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(r,&value));
19173c0623SSatish Balay     rand1[i] = (int)(value* 144327);
205f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(r,&value));
21173c0623SSatish Balay     rand2[i] = (int)(value* 144327);
22173c0623SSatish Balay   }
23173c0623SSatish Balay 
2498d6af09SSatish Balay   /* Take care of paging effects */
255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(100,&dummy));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(dummy));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTime(&x));
2898d6af09SSatish Balay 
29173c0623SSatish Balay   /* Do all mallocs */
30173c0623SSatish Balay   for (i=0; i< 1000; i++) {
315f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(rand1[i],&arr[i]));
32173c0623SSatish Balay   }
33173c0623SSatish Balay 
345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTime(&x));
35173c0623SSatish Balay 
36173c0623SSatish Balay   /* Do some frees */
37173c0623SSatish Balay   for (i=0; i< 1000; i+=2) {
385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(arr[i]));
39173c0623SSatish Balay   }
40173c0623SSatish Balay 
41173c0623SSatish Balay   /* Do some mallocs */
42173c0623SSatish Balay   for (i=0; i< 1000; i+=2) {
435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(rand2[i],&arr[i]));
44173c0623SSatish Balay   }
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTime(&y));
46173c0623SSatish Balay 
47173c0623SSatish Balay   for (i=0; i< 1000; i++) {
485f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(arr[i]));
49173c0623SSatish Balay   }
50173c0623SSatish Balay 
51b4d8b9abSSatish Balay   fprintf(stdout,"%-15s : %e sec, with options : ","PetscMalloc",(y-x)/500.0);
525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,"-malloc",&flg));
538caf3d72SBarry Smith   if (flg) fprintf(stdout,"-malloc ");
54b4d8b9abSSatish Balay   fprintf(stdout,"\n");
55173c0623SSatish Balay 
565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&r));
57*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
58*b122ec5aSJacob Faibussowitsch   return 0;
59173c0623SSatish Balay }
60