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