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