Actual source code: ex52.c
petsc-3.13.6 2020-09-29
1: static char help[] = "A benchmark for testing PetscSortInt() and PetscSortIntWithArrayPair()\n\
2: The array is filled with random numbers, but one can control average duplicates for each unique integer with the -d option.\n\
3: Usage:\n\
4: mpirun -n 1 ./ex52 -n <length of the array to sort>, default=100 \n\
5: -r <repeat times for each sort>, default=10 \n\
6: -d <average duplicates for each unique integer>, default=1, i.e., no duplicates \n\n";
8: #include <petscsys.h>
9: #include <petsctime.h>
10: int main(int argc,char **argv)
11: {
13: PetscInt i,l,n=100,r=10,d=1;
14: PetscInt *X,*Y,*Z;
15: PetscReal val;
16: PetscRandom rdm;
17: PetscLogDouble time;
18: PetscMPIInt size;
20: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
21: MPI_Comm_size(PETSC_COMM_WORLD,&size);
22: if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"This is a uniprocessor example only!");
24: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
25: PetscOptionsGetInt(NULL,NULL,"-r",&r,NULL);
26: PetscOptionsGetInt(NULL,NULL,"-d",&d,NULL);
27: if (n<1 || r<1 || d<1 || d>n) SETERRQ3(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Wrong input n=%D,r=%D,d=%d. They must be >=1 and n>=d\n",n,r,d);
29: PetscCalloc3(n,&X,n,&Y,n,&Z);
30: PetscRandomCreate(PETSC_COMM_SELF,&rdm);
31: PetscRandomSetFromOptions(rdm);
33: time = 0.0;
34: for (l=0; l<r; l++) { /* r loops */
35: for (i=0; i<n; i++) { /* Init X[] */
36: PetscRandomGetValueReal(rdm,&val);
37: X[i] = val*PETSC_MAX_INT;
38: if (d > 1) X[i] = X[i] % (n/d);
39: }
41: PetscTimeSubtract(&time);
42: PetscSortInt(n,X);
43: PetscTimeAdd(&time);
45: for (i=0; i<n-1; i++) {if (X[i] > X[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PetscSortInt() produced wrong results!");}
46: }
47: PetscPrintf(PETSC_COMM_SELF,"PetscSortInt() with %D integers, %D duplicate(s) per unique value took %g seconds\n",n,d,time/r);
49: time = 0.0;
50: for (l=0; l<r; l++) { /* r loops */
51: for (i=0; i<n; i++) { /* Init X[] */
52: PetscRandomGetValueReal(rdm,&val);
53: X[i] = val*PETSC_MAX_INT;
54: if (d > 1) X[i] = X[i] % (n/d);
55: }
57: PetscTimeSubtract(&time);
58: PetscSortIntWithArrayPair(n,X,Y,Z);
59: PetscTimeAdd(&time);
61: for (i=0; i<n-1; i++) {if (X[i] > X[i+1]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"PetscSortInt() produced wrong results!");}
62: }
63: PetscPrintf(PETSC_COMM_SELF,"PetscSortIntWithArrayPair() with %D integers, %D duplicate(s) per unique value took %g seconds\n",n,d,time/r);
64: PetscPrintf(PETSC_COMM_SELF,"SUCCEEDED\n");
66: PetscRandomDestroy(&rdm);
67: PetscFree3(X,Y,Z);
68: PetscFinalize();
69: return ierr;
70: }
72: /*TEST
74: test:
75: args: -n 1000 -r 1 -d 1
76: # Do not need to output timing results for test
77: filter: grep -v "per unique value took"
79: TEST*/