xref: /petsc/src/sys/tests/ex12.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[] = "Tests timing PetscSortInt().\n\n";
3 
4 #include <petscsys.h>
5 
6 int main(int argc,char **argv)
7 {
8   PetscInt       i,n = 1000,*values;
9 #if defined(PETSC_USE_LOG)
10   PetscLogEvent  event;
11 #endif
12   PetscRandom    rand;
13   PetscReal      value;
14   PetscErrorCode ierr;
15   PetscBool      values_view=PETSC_FALSE;
16   PetscMPIInt    rank;
17 
18   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
19   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
20   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
21   CHKERRQ(PetscOptionsGetBool(NULL,0,"-values_view",&values_view,NULL));
22 
23   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&rand));
24   CHKERRQ(PetscRandomSetFromOptions(rand));
25 
26   CHKERRQ(PetscMalloc1(n,&values));
27   for (i=0; i<n; i++) {
28     CHKERRQ(PetscRandomGetValueReal(rand,&value));
29     values[i] = (PetscInt)(n*value + 2.0);
30   }
31   CHKERRQ(PetscSortInt(n,values));
32 
33   CHKERRQ(PetscLogEventRegister("Sort",0,&event));
34   CHKERRQ(PetscLogEventBegin(event,0,0,0,0));
35 
36   for (i=0; i<n; i++) {
37     CHKERRQ(PetscRandomGetValueReal(rand,&value));
38     values[i] = (PetscInt)(n*value + 2.0);
39   }
40   CHKERRQ(PetscSortInt(n,values));
41   CHKERRQ(PetscLogEventEnd(event,0,0,0,0));
42 
43   for (i=1; i<n; i++) {
44     PetscCheckFalse(values[i] < values[i-1],PETSC_COMM_SELF,PETSC_ERR_PLIB,"Values not sorted");
45     if (values_view && rank == 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT " %" PetscInt_FMT "\n",i,values[i]));
46   }
47   CHKERRQ(PetscFree(values));
48   CHKERRQ(PetscRandomDestroy(&rand));
49 
50   ierr = PetscFinalize();
51   return ierr;
52 }
53 
54 /*TEST
55 
56    test:
57       args: -values_view
58 
59 TEST*/
60