xref: /petsc/src/benchmarks/PetscMalloc.c (revision b4d8b9abfc967e8f3c9331284eb41aa08b1075bc)
1*b4d8b9abSSatish Balay /*$Id: PetscMalloc.c,v 1.27 2001/08/29 20:32:42 balay Exp balay $*/
2d8e9fea7SSatish Balay 
3173c0623SSatish Balay #include "petsc.h"
4e090d566SSatish Balay #include "petscsys.h"
5173c0623SSatish Balay 
64a2ae208SSatish Balay #undef __FUNCT__
74a2ae208SSatish Balay #define __FUNCT__ "main"
8173c0623SSatish Balay int main(int argc,char **argv)
9173c0623SSatish Balay {
10b0a32e0cSBarry Smith   PetscLogDouble  x,y;
1147794344SBarry Smith   double      value;
1298d6af09SSatish Balay   void        *arr[1000],*dummy;
13f1af5d2fSBarry Smith   int         ierr,i,rand1[1000],rand2[1000];
1477c4ece6SBarry Smith   PetscRandom r;
15f1af5d2fSBarry Smith   PetscTruth  flg;
16173c0623SSatish Balay 
1777c4ece6SBarry Smith   PetscInitialize(&argc,&argv,0,0);
18173c0623SSatish Balay 
19029af93fSBarry Smith   ierr = PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&r);CHKERRQ(ierr);
20173c0623SSatish Balay   for (i=0; i<1000; i++) {
2177c4ece6SBarry Smith     ierr    = PetscRandomGetValue(r,&value);CHKERRQ(ierr);
22173c0623SSatish Balay     rand1[i] = (int)(value* 144327);
2377c4ece6SBarry Smith     ierr    = PetscRandomGetValue(r,&value);CHKERRQ(ierr);
24173c0623SSatish Balay     rand2[i] = (int)(value* 144327);
25173c0623SSatish Balay   }
26173c0623SSatish Balay 
2798d6af09SSatish Balay   /* Take care of paging effects */
28ac355199SBarry Smith   ierr = PetscMalloc(100,&dummy);CHKERRQ(ierr);
29ac355199SBarry Smith   ierr = PetscFree(dummy);CHKERRQ(ierr);
30ac355199SBarry Smith   ierr = PetscGetTime(&x);CHKERRQ(ierr);
3198d6af09SSatish Balay 
32173c0623SSatish Balay   /* Do all mallocs */
33173c0623SSatish Balay   for (i=0 ; i< 1000; i++) {
34ac355199SBarry Smith     ierr = PetscMalloc(rand1[i],& arr[i]);CHKERRQ(ierr);
35173c0623SSatish Balay   }
36173c0623SSatish Balay 
37ac355199SBarry Smith   ierr = PetscGetTime(&x);CHKERRQ(ierr);
38173c0623SSatish Balay 
39173c0623SSatish Balay   /* Do some frees */
40173c0623SSatish Balay   for (i=0; i< 1000; i+=2) {
41ac355199SBarry Smith     ierr = PetscFree(arr[i]);CHKERRQ(ierr);
42173c0623SSatish Balay   }
43173c0623SSatish Balay 
44173c0623SSatish Balay   /* Do some mallocs */
45173c0623SSatish Balay   for (i=0; i< 1000; i+=2) {
46ac355199SBarry Smith     ierr = PetscMalloc(rand2[i],&arr[i]);CHKERRQ(ierr);
47173c0623SSatish Balay  }
48ac355199SBarry Smith   ierr = PetscGetTime(&y);CHKERRQ(ierr);
49173c0623SSatish Balay 
50173c0623SSatish Balay   for (i=0; i< 1000; i++) {
51ac355199SBarry Smith     ierr = PetscFree(arr[i]);CHKERRQ(ierr);
52173c0623SSatish Balay   }
53173c0623SSatish Balay 
54*b4d8b9abSSatish Balay   fprintf(stdout,"%-15s : %e sec, with options : ","PetscMalloc",(y-x)/500.0);
55*b4d8b9abSSatish Balay   if(PetscOptionsHasName(PETSC_NULL,"-trmalloc",&flg),flg) fprintf(stdout,"-trmalloc ");
56*b4d8b9abSSatish Balay   fprintf(stdout,"\n");
57173c0623SSatish Balay 
58ac355199SBarry Smith   ierr = PetscRandomDestroy(r);CHKERRQ(ierr);
592b7fea2aSSatish Balay   ierr = PetscFinalize();CHKERRQ(ierr);
603a40ed3dSBarry Smith   PetscFunctionReturn(0);
61173c0623SSatish Balay }
62