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