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