xref: /petsc/src/benchmarks/PetscMemcmp.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1d8e9fea7SSatish Balay 
2c6db04a5SJed Brown #include <petscsys.h>
38563dfccSBarry Smith #include <petsctime.h>
4b0878937SSatish Balay 
5b0878937SSatish Balay int main(int argc,char **argv)
6b0878937SSatish Balay {
7b0a32e0cSBarry Smith   PetscLogDouble x,y,z;
8ea709b57SSatish Balay   PetscScalar    A[10000],B[10000];
9a438ae71SBarry Smith   PetscInt       i;
10ace3abfcSBarry Smith   PetscBool      flg;
11b0878937SSatish Balay 
12*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,0,0));
13b0878937SSatish Balay   for (i=0; i<10000; i++) {
14b0878937SSatish Balay     A[i] = i%61897;
15b0878937SSatish Balay     B[i] = i%61897;
16b0878937SSatish Balay   }
17b0878937SSatish Balay   /* To take care of paging effects */
184115e8b9SSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*0,&flg);
195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTime(&x));
20b0878937SSatish Balay 
215f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTime(&x));
224115e8b9SSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*10000,&flg);
234115e8b9SSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*10000,&flg);
244115e8b9SSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*10000,&flg);
254115e8b9SSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*10000,&flg);
264115e8b9SSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*10000,&flg);
274115e8b9SSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*10000,&flg);
284115e8b9SSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*10000,&flg);
294115e8b9SSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*10000,&flg);
304115e8b9SSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*10000,&flg);
314115e8b9SSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*10000,&flg);
325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTime(&y));
339ae0b57aSSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*0,&flg);
349ae0b57aSSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*0,&flg);
359ae0b57aSSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*0,&flg);
369ae0b57aSSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*0,&flg);
379ae0b57aSSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*0,&flg);
389ae0b57aSSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*0,&flg);
399ae0b57aSSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*0,&flg);
409ae0b57aSSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*0,&flg);
419ae0b57aSSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*0,&flg);
429ae0b57aSSatish Balay   PetscMemcmp(A,B,sizeof(PetscScalar)*0,&flg);
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscTime(&z));
44b0878937SSatish Balay 
4535d8aa7fSBarry Smith   fprintf(stdout,"%s : \n","PetscMemcmp");
46b4d8b9abSSatish Balay   fprintf(stdout,"    %-15s : %e sec\n","Latency",(z-y)/10.0);
47b4d8b9abSSatish Balay   fprintf(stdout,"    %-15s : %e sec\n","Per PetscScalar",(2*y-x-z)/100000);
48b0878937SSatish Balay 
49*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
50*b122ec5aSJacob Faibussowitsch   return 0;
51b0878937SSatish Balay }
52