xref: /petsc/src/sys/tests/ex35.c (revision 9566063d113dddea24716c546802770db7481bc0)
1c4762a1bSJed Brown static char help[] = "Tests PetscSortReal(), PetscSortRealWithArrayInt(), PetscFindReal\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscsys.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc,char **argv)
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   PetscInt       i,loc;
8c4762a1bSJed Brown   PetscReal      val;
9c4762a1bSJed Brown   PetscReal      x[] = {39, 9, 19, -39, 29, 309, 209, 390, 12, 11};
10c4762a1bSJed Brown   PetscReal      x2[] = {39, 9, 19, -39, 29, 309, 209, 390, 12, 11};
11c4762a1bSJed Brown   PetscReal      x3[] = {39, 9, 19, -39, 29, 309, 209, 390, 12, 11};
12c4762a1bSJed Brown   PetscInt       index2[] = {1, -1, 4, 12, 13, 14, 0, 7, 9, 11};
13c4762a1bSJed Brown   PetscInt       index3[] = {1, -1, 4, 12, 13, 14, 0, 7, 9, 11};
14c4762a1bSJed Brown 
15*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
16*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"1st test\n"));
17*9566063dSJacob Faibussowitsch   for (i=0; i<5; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF," %g\n",(double)x[i]));
18*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"---------------\n"));
19*9566063dSJacob Faibussowitsch   PetscCall(PetscSortReal(5,x));
20*9566063dSJacob Faibussowitsch   for (i=0; i<5; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF," %g\n",(double)x[i]));
21c4762a1bSJed Brown 
22*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n2nd test\n"));
23*9566063dSJacob Faibussowitsch   for (i=0; i<10; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF," %g\n",(double)x[i]));
24*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"---------------\n"));
25*9566063dSJacob Faibussowitsch   PetscCall(PetscSortReal(10,x));
26*9566063dSJacob Faibussowitsch   for (i=0; i<10; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF," %g\n",(double)x[i]));
27c4762a1bSJed Brown 
28*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n3rd test\n"));
29*9566063dSJacob Faibussowitsch   for (i=0; i<5; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF," %2" PetscInt_FMT "     %g\n",index2[i], (double)x2[i]));
30*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"---------------\n"));
31*9566063dSJacob Faibussowitsch   PetscCall(PetscSortRealWithArrayInt(5, x2, index2));
32*9566063dSJacob Faibussowitsch   for (i=0; i<5; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF," %2" PetscInt_FMT "     %g\n",index2[i], (double)x2[i]));
33c4762a1bSJed Brown 
34*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n4th test\n"));
35*9566063dSJacob Faibussowitsch   for (i=0; i<10; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF," %2" PetscInt_FMT "     %g\n",index3[i], (double)x3[i]));
36*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"---------------\n"));
37*9566063dSJacob Faibussowitsch   PetscCall(PetscSortRealWithArrayInt(10, x3, index3));
38*9566063dSJacob Faibussowitsch   for (i=0; i<10; i++) PetscCall(PetscPrintf(PETSC_COMM_SELF," %2" PetscInt_FMT "     %g\n",index3[i], (double)x3[i]));
39c4762a1bSJed Brown 
40*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n5th test\n"));
41c4762a1bSJed Brown   val  = 44;
42*9566063dSJacob Faibussowitsch   PetscCall(PetscFindReal(val,10,x3,PETSC_SMALL,&loc));
43*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF," %g in array: loc %" PetscInt_FMT "\n",(double)val,loc));
44c4762a1bSJed Brown   val  = 309.2;
45*9566063dSJacob Faibussowitsch   PetscCall(PetscFindReal(val,10,x3,PETSC_SMALL,&loc));
46*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF," %g in array: loc %" PetscInt_FMT "\n",(double)val,loc));
47c4762a1bSJed Brown   val  = 309.2;
48*9566063dSJacob Faibussowitsch   PetscCall(PetscFindReal(val,10,x3,0.21,&loc));
49*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF," %g in array: loc %" PetscInt_FMT "\n",(double)val,loc));
50c4762a1bSJed Brown 
51*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
52b122ec5aSJacob Faibussowitsch   return 0;
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown /*TEST
56c4762a1bSJed Brown 
57c4762a1bSJed Brown    test:
58c4762a1bSJed Brown 
59c4762a1bSJed Brown TEST*/
60