Actual source code: PetscMemcpy.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2:  #include <petscsys.h>
  3:  #include <petsctime.h>

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

 12:   PetscInitialize(&argc,&argv,0,0);if (ierr) return ierr;
 13:   PetscMalloc1(8000000,&A);
 14:   PetscMalloc1(8000000,&B);

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

 24:   PetscTime(&x);
 25:   /*
 26:   PetscMemcpy(A,B,sizeof(PetscScalar)*8000000);
 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:   { int j;
 37:     for (j = 0; j<10; j++)
 38:       for (i=0; i<8000000; i++) B[i] = A[i];}

 40:   PetscTime(&y);
 41:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 42:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 43:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 44:   PetscMemcpy(A,B,sizeof(PetscScalar)*0);
 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:   PetscTime(&z);

 53:   fprintf(stdout,"%s : \n","PetscMemcpy");
 54:   fprintf(stdout,"    %-15s : %e MB/s\n","Bandwidth",10.0*8*8/(y-x));
 55:   fprintf(stdout,"    %-15s : %e sec\n","Latency",(z-y)/10.0);
 56:   fprintf(stdout,"    %-15s : %e sec\n","Per PetscScalar",(2*y-x-z)/8000000.0);

 58:   PetscFree(A);
 59:   PetscFree(B);
 60:   PetscFinalize();
 61:   return ierr;
 62: }