1 2 static char help[] = "Tests PC and KSP on a tridiagonal matrix. Note that most\n\ 3 users should employ the KSP interface instead of using PC directly.\n\n"; 4 5 6 #include <petscksp.h> 7 8 int main(int argc,char **args) 9 { 10 Mat mat; /* matrix */ 11 Vec b,ustar,u; /* vectors (RHS, exact solution, approx solution) */ 12 PC pc; /* PC context */ 13 KSP ksp; /* KSP context */ 14 PetscErrorCode ierr; 15 PetscInt n = 10,i,its,col[3]; 16 PetscScalar value[3]; 17 PCType pcname; 18 KSPType kspname; 19 PetscReal norm,tol=1000.*PETSC_MACHINE_EPSILON; 20 21 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 22 /* Create and initialize vectors */ 23 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&b);CHKERRQ(ierr); 24 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&ustar);CHKERRQ(ierr); 25 ierr = VecCreateSeq(PETSC_COMM_SELF,n,&u);CHKERRQ(ierr); 26 ierr = VecSet(ustar,1.0);CHKERRQ(ierr); 27 ierr = VecSet(u,0.0);CHKERRQ(ierr); 28 29 /* Create and assemble matrix */ 30 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&mat);CHKERRQ(ierr); 31 value[0] = -1.0; value[1] = 2.0; value[2] = -1.0; 32 for (i=1; i<n-1; i++) { 33 col[0] = i-1; col[1] = i; col[2] = i+1; 34 ierr = MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 35 } 36 i = n - 1; col[0] = n - 2; col[1] = n - 1; 37 ierr = MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 38 i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0; 39 ierr = MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 40 ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 41 ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 42 43 /* Compute right-hand-side vector */ 44 ierr = MatMult(mat,ustar,b);CHKERRQ(ierr); 45 46 /* Create PC context and set up data structures */ 47 ierr = PCCreate(PETSC_COMM_WORLD,&pc);CHKERRQ(ierr); 48 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 49 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); 50 ierr = PCSetOperators(pc,mat,mat);CHKERRQ(ierr); 51 ierr = PCSetUp(pc);CHKERRQ(ierr); 52 53 /* Create KSP context and set up data structures */ 54 ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr); 55 ierr = KSPSetType(ksp,KSPRICHARDSON);CHKERRQ(ierr); 56 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 57 ierr = PCSetOperators(pc,mat,mat);CHKERRQ(ierr); 58 ierr = KSPSetPC(ksp,pc);CHKERRQ(ierr); 59 ierr = KSPSetUp(ksp);CHKERRQ(ierr); 60 61 /* Solve the problem */ 62 ierr = KSPGetType(ksp,&kspname);CHKERRQ(ierr); 63 ierr = PCGetType(pc,&pcname);CHKERRQ(ierr); 64 ierr = PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioning\n",kspname,pcname);CHKERRQ(ierr); 65 ierr = KSPSolve(ksp,b,u);CHKERRQ(ierr); 66 ierr = VecAXPY(u,-1.0,ustar);CHKERRQ(ierr); 67 ierr = VecNorm(u,NORM_2,&norm);CHKERRQ(ierr); 68 ierr = KSPGetIterationNumber(ksp,&its);CHKERRQ(ierr); 69 if (norm > tol) { 70 ierr = PetscPrintf(PETSC_COMM_SELF,"2 norm of error %g Number of iterations %D\n",(double)norm,its);CHKERRQ(ierr); 71 } 72 73 /* Free data structures */ 74 ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 75 ierr = VecDestroy(&u);CHKERRQ(ierr); 76 ierr = VecDestroy(&ustar);CHKERRQ(ierr); 77 ierr = VecDestroy(&b);CHKERRQ(ierr); 78 ierr = MatDestroy(&mat);CHKERRQ(ierr); 79 ierr = PCDestroy(&pc);CHKERRQ(ierr); 80 81 ierr = PetscFinalize(); 82 return ierr; 83 } 84 85 86 87 88 89 /*TEST 90 91 test: 92 args: -ksp_type cg -ksp_monitor_short 93 94 TEST*/ 95