xref: /petsc/src/sys/tests/ex12.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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   PetscBool      values_view=PETSC_FALSE;
15c4762a1bSJed Brown   PetscMPIInt    rank;
16c4762a1bSJed Brown 
179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,0,"-values_view",&values_view,NULL));
21c4762a1bSJed Brown 
229566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rand));
239566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(rand));
24c4762a1bSJed Brown 
259566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n,&values));
26c4762a1bSJed Brown   for (i=0; i<n; i++) {
279566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(rand,&value));
28c4762a1bSJed Brown     values[i] = (PetscInt)(n*value + 2.0);
29c4762a1bSJed Brown   }
309566063dSJacob Faibussowitsch   PetscCall(PetscSortInt(n,values));
31c4762a1bSJed Brown 
329566063dSJacob Faibussowitsch   PetscCall(PetscLogEventRegister("Sort",0,&event));
339566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(event,0,0,0,0));
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   for (i=0; i<n; i++) {
369566063dSJacob Faibussowitsch     PetscCall(PetscRandomGetValueReal(rand,&value));
37c4762a1bSJed Brown     values[i] = (PetscInt)(n*value + 2.0);
38c4762a1bSJed Brown   }
399566063dSJacob Faibussowitsch   PetscCall(PetscSortInt(n,values));
409566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(event,0,0,0,0));
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   for (i=1; i<n; i++) {
43*08401ef6SPierre Jolivet     PetscCheck(values[i] >= values[i-1],PETSC_COMM_SELF,PETSC_ERR_PLIB,"Values not sorted");
449566063dSJacob Faibussowitsch     if (values_view && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF,"%" PetscInt_FMT " %" PetscInt_FMT "\n",i,values[i]));
45c4762a1bSJed Brown   }
469566063dSJacob Faibussowitsch   PetscCall(PetscFree(values));
479566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&rand));
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
50b122ec5aSJacob Faibussowitsch   return 0;
51c4762a1bSJed Brown }
52c4762a1bSJed Brown 
53c4762a1bSJed Brown /*TEST
54c4762a1bSJed Brown 
55c4762a1bSJed Brown    test:
56c4762a1bSJed Brown       args: -values_view
57c4762a1bSJed Brown 
58c4762a1bSJed Brown TEST*/
59