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