Actual source code: PetscMemcpy.c

petsc-3.3-p7 2013-05-11
  2: #include <petscsys.h>

  6: int main(int argc,char **argv)
  7: {
  8:   PetscLogDouble  x,y,z;
  9:   int             i,ierr;
 10:   PetscScalar     *A,*B;

 12:   PetscInitialize(&argc,&argv,0,0);

 14:   PetscMalloc(8000000*sizeof(PetscScalar),&A);
 15:   PetscMalloc(8000000*sizeof(PetscScalar),&B);

 17:   for (i=0; i<8000000; i++) {
 18:     A[i] = i%61897;
 19:     B[i] = i%61897;
 20:   }
 21:   /* To take care of paging effects */
 22:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 23:   PetscGetTime(&x);

 25:   PetscGetTime(&x);
 26:   /*
 27:   PetscMemcpy(A,B,sizeof(PetscScalar)*8000000);
 28:   PetscMemcpy(A,B,sizeof(PetscScalar)*8000000);
 29:   PetscMemcpy(A,B,sizeof(PetscScalar)*8000000);
 30:   PetscMemcpy(A,B,sizeof(PetscScalar)*8000000);
 31:   PetscMemcpy(A,B,sizeof(PetscScalar)*8000000);
 32:   PetscMemcpy(A,B,sizeof(PetscScalar)*8000000);
 33:   PetscMemcpy(A,B,sizeof(PetscScalar)*8000000);
 34:   PetscMemcpy(A,B,sizeof(PetscScalar)*8000000);
 35:   PetscMemcpy(A,B,sizeof(PetscScalar)*8000000);
 36:   PetscMemcpy(A,B,sizeof(PetscScalar)*8000000); */
 37:   { int j;
 38:   for (j = 0; j<10; j++) {
 39:     for (i=0; i<8000000; i++) {
 40:       B[i] = A[i];
 41:     }
 42:   }}

 44:   PetscGetTime(&y);
 45:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 46:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 47:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 48:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 49:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 50:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 51:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 52:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 53:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 54:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 55:   PetscGetTime(&z);

 57:   fprintf(stdout,"%s : \n","PetscMemcpy");
 58:   fprintf(stdout,"    %-15s : %e MB/s\n","Bandwidth",10.0*8*8/(y-x));
 59:   fprintf(stdout,"    %-15s : %e sec\n","Latency",(z-y)/10.0);
 60:   fprintf(stdout,"    %-15s : %e sec\n","Per PetscScalar",(2*y-x-z)/8000000.0);

 62:   PetscFree(A);
 63:   PetscFree(B);
 64:   PetscFinalize();
 65:   return(0);
 66: }