1c4762a1bSJed Brown static char help[] = "Tests timing PetscSortInt().\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscsys.h>
4c4762a1bSJed Brown
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown PetscInt i, n = 1000, *values;
8956f8c0dSBarry Smith PetscLogEvent event;
9c4762a1bSJed Brown PetscRandom rand;
10c4762a1bSJed Brown PetscReal value;
11c4762a1bSJed Brown PetscBool values_view = PETSC_FALSE;
12c4762a1bSJed Brown PetscMPIInt rank;
13c4762a1bSJed Brown
14327415f7SBarry Smith PetscFunctionBeginUser;
15*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-n", &n, NULL));
189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, 0, "-values_view", &values_view, NULL));
19c4762a1bSJed Brown
209566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rand));
219566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand));
22c4762a1bSJed Brown
239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &values));
24c4762a1bSJed Brown for (i = 0; i < n; i++) {
259566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &value));
26c4762a1bSJed Brown values[i] = (PetscInt)(n * value + 2.0);
27c4762a1bSJed Brown }
289566063dSJacob Faibussowitsch PetscCall(PetscSortInt(n, values));
29c4762a1bSJed Brown
309566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("Sort", 0, &event));
319566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(event, 0, 0, 0, 0));
32c4762a1bSJed Brown
33c4762a1bSJed Brown for (i = 0; i < n; i++) {
349566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(rand, &value));
35c4762a1bSJed Brown values[i] = (PetscInt)(n * value + 2.0);
36c4762a1bSJed Brown }
379566063dSJacob Faibussowitsch PetscCall(PetscSortInt(n, values));
389566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(event, 0, 0, 0, 0));
39c4762a1bSJed Brown
40c4762a1bSJed Brown for (i = 1; i < n; i++) {
4108401ef6SPierre Jolivet PetscCheck(values[i] >= values[i - 1], PETSC_COMM_SELF, PETSC_ERR_PLIB, "Values not sorted");
429566063dSJacob Faibussowitsch if (values_view && rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT " %" PetscInt_FMT "\n", i, values[i]));
43c4762a1bSJed Brown }
449566063dSJacob Faibussowitsch PetscCall(PetscFree(values));
459566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand));
46c4762a1bSJed Brown
479566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
48b122ec5aSJacob Faibussowitsch return 0;
49c4762a1bSJed Brown }
50c4762a1bSJed Brown
51c4762a1bSJed Brown /*TEST
52c4762a1bSJed Brown
53c4762a1bSJed Brown test:
54c4762a1bSJed Brown args: -values_view
55c4762a1bSJed Brown
56c4762a1bSJed Brown TEST*/
57