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 PetscErrorCode ierr; 15 16 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 17 ierr = PetscPrintf(PETSC_COMM_SELF,"1st test\n");CHKERRQ(ierr); 18 for (i=0; i<5; i++) {ierr = PetscPrintf(PETSC_COMM_SELF," %g\n",(double)x[i]);CHKERRQ(ierr);} 19 ierr = PetscPrintf(PETSC_COMM_SELF,"---------------\n");CHKERRQ(ierr); 20 ierr = PetscSortReal(5,x);CHKERRQ(ierr); 21 for (i=0; i<5; i++) {ierr = PetscPrintf(PETSC_COMM_SELF," %g\n",(double)x[i]);CHKERRQ(ierr);} 22 23 ierr = PetscPrintf(PETSC_COMM_SELF,"\n2nd test\n");CHKERRQ(ierr); 24 for (i=0; i<10; i++) {ierr = PetscPrintf(PETSC_COMM_SELF," %g\n",(double)x[i]);CHKERRQ(ierr);} 25 ierr = PetscPrintf(PETSC_COMM_SELF,"---------------\n");CHKERRQ(ierr); 26 ierr = PetscSortReal(10,x);CHKERRQ(ierr); 27 for (i=0; i<10; i++) {ierr = PetscPrintf(PETSC_COMM_SELF," %g\n",(double)x[i]);CHKERRQ(ierr);} 28 29 ierr = PetscPrintf(PETSC_COMM_SELF,"\n3rd test\n");CHKERRQ(ierr); 30 for (i=0; i<5; i++) {ierr = PetscPrintf(PETSC_COMM_SELF," %2D %g\n",index2[i], (double)x2[i]);CHKERRQ(ierr);} 31 ierr = PetscPrintf(PETSC_COMM_SELF,"---------------\n");CHKERRQ(ierr); 32 ierr = PetscSortRealWithArrayInt(5, x2, index2);CHKERRQ(ierr); 33 for (i=0; i<5; i++) {ierr = PetscPrintf(PETSC_COMM_SELF," %2D %g\n",index2[i], (double)x2[i]);CHKERRQ(ierr);} 34 35 ierr = PetscPrintf(PETSC_COMM_SELF,"\n4th test\n");CHKERRQ(ierr); 36 for (i=0; i<10; i++) {ierr = PetscPrintf(PETSC_COMM_SELF," %2D %g\n",index3[i], (double)x3[i]);CHKERRQ(ierr);} 37 ierr = PetscPrintf(PETSC_COMM_SELF,"---------------\n");CHKERRQ(ierr); 38 ierr = PetscSortRealWithArrayInt(10, x3, index3);CHKERRQ(ierr); 39 for (i=0; i<10; i++) {ierr = PetscPrintf(PETSC_COMM_SELF," %2D %g\n",index3[i], (double)x3[i]);CHKERRQ(ierr);} 40 41 ierr = PetscPrintf(PETSC_COMM_SELF,"\n5th test\n");CHKERRQ(ierr); 42 val = 44; 43 ierr = PetscFindReal(val,10,x3,PETSC_SMALL,&loc);CHKERRQ(ierr); 44 ierr = PetscPrintf(PETSC_COMM_SELF," %g in array: loc %D\n",(double)val,loc);CHKERRQ(ierr); 45 val = 309.2; 46 ierr = PetscFindReal(val,10,x3,PETSC_SMALL,&loc);CHKERRQ(ierr); 47 ierr = PetscPrintf(PETSC_COMM_SELF," %g in array: loc %D\n",(double)val,loc);CHKERRQ(ierr); 48 val = 309.2; 49 ierr = PetscFindReal(val,10,x3,0.21,&loc);CHKERRQ(ierr); 50 ierr = PetscPrintf(PETSC_COMM_SELF," %g in array: loc %D\n",(double)val,loc);CHKERRQ(ierr); 51 52 ierr = PetscFinalize(); 53 return ierr; 54 } 55 56 /*TEST 57 58 test: 59 60 61 TEST*/ 62