xref: /petsc/src/benchmarks/PetscMalloc.c (revision abdd934af67c1cc360c4d0cfcf2de16e188d627e)
1c6db04a5SJed Brown #include <petscsys.h>
28563dfccSBarry Smith #include <petsctime.h>
3173c0623SSatish Balay 
main(int argc,char ** argv)4173c0623SSatish Balay int main(int argc, char **argv)
5173c0623SSatish Balay {
6b0a32e0cSBarry Smith   PetscLogDouble x, y;
747794344SBarry Smith   double         value;
898d6af09SSatish Balay   void          *arr[1000], *dummy;
9a438ae71SBarry Smith   int            i, rand1[1000], rand2[1000];
1077c4ece6SBarry Smith   PetscRandom    r;
11ace3abfcSBarry Smith   PetscBool      flg;
12173c0623SSatish Balay 
139566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, 0, 0));
149566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r));
159566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(r));
16173c0623SSatish Balay   for (i = 0; i < 1000; i++) {
179566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(r, &value));
18173c0623SSatish Balay     rand1[i] = (int)(value * 144327);
199566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValue(r, &value));
20173c0623SSatish Balay     rand2[i] = (int)(value * 144327);
21173c0623SSatish Balay   }
22173c0623SSatish Balay 
2398d6af09SSatish Balay   /* Take care of paging effects */
249566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(100, &dummy));
259566063dSJacob Faibussowitsch   PetscCall(PetscFree(dummy));
269566063dSJacob Faibussowitsch   PetscCall(PetscTime(&x));
2798d6af09SSatish Balay 
28173c0623SSatish Balay   /* Do all mallocs */
29*3a7d0413SPierre Jolivet   for (i = 0; i < 1000; i++) PetscCall(PetscMalloc1(rand1[i], &arr[i]));
30173c0623SSatish Balay 
319566063dSJacob Faibussowitsch   PetscCall(PetscTime(&x));
32173c0623SSatish Balay 
33173c0623SSatish Balay   /* Do some frees */
34*3a7d0413SPierre Jolivet   for (i = 0; i < 1000; i += 2) PetscCall(PetscFree(arr[i]));
35173c0623SSatish Balay 
36173c0623SSatish Balay   /* Do some mallocs */
37*3a7d0413SPierre Jolivet   for (i = 0; i < 1000; i += 2) PetscCall(PetscMalloc1(rand2[i], &arr[i]));
389566063dSJacob Faibussowitsch   PetscCall(PetscTime(&y));
39173c0623SSatish Balay 
40*3a7d0413SPierre Jolivet   for (i = 0; i < 1000; i++) PetscCall(PetscFree(arr[i]));
41173c0623SSatish Balay 
42b4d8b9abSSatish Balay   fprintf(stdout, "%-15s : %e sec, with options : ", "PetscMalloc", (y - x) / 500.0);
439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, "-malloc", &flg));
448caf3d72SBarry Smith   if (flg) fprintf(stdout, "-malloc ");
45b4d8b9abSSatish Balay   fprintf(stdout, "\n");
46173c0623SSatish Balay 
479566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&r));
489566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
49b122ec5aSJacob Faibussowitsch   return 0;
50173c0623SSatish Balay }
51