Actual source code: ex2.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  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";

  5: #include <petscksp.h>

  9: int main(int argc,char **args)
 10: {
 11:   Mat            mat;          /* matrix */
 12:   Vec            b,ustar,u;  /* vectors (RHS, exact solution, approx solution) */
 13:   PC             pc;           /* PC context */
 14:   KSP            ksp;          /* KSP context */
 16:   PetscInt       n = 10,i,its,col[3];
 17:   PetscScalar    value[3];
 18:   PCType         pcname;
 19:   KSPType        kspname;
 20:   PetscReal      norm,tol=1.e-14;

 22:   PetscInitialize(&argc,&args,(char*)0,help);

 24:   /* Create and initialize vectors */
 25:   VecCreateSeq(PETSC_COMM_SELF,n,&b);
 26:   VecCreateSeq(PETSC_COMM_SELF,n,&ustar);
 27:   VecCreateSeq(PETSC_COMM_SELF,n,&u);
 28:   VecSet(ustar,1.0);
 29:   VecSet(u,0.0);

 31:   /* Create and assemble matrix */
 32:   MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&mat);
 33:   value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 34:   for (i=1; i<n-1; i++) {
 35:     col[0] = i-1; col[1] = i; col[2] = i+1;
 36:     MatSetValues(mat,1,&i,3,col,value,INSERT_VALUES);
 37:   }
 38:   i    = n - 1; col[0] = n - 2; col[1] = n - 1;
 39:   MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);
 40:   i    = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
 41:   MatSetValues(mat,1,&i,2,col,value,INSERT_VALUES);
 42:   MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);
 43:   MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);

 45:   /* Compute right-hand-side vector */
 46:   MatMult(mat,ustar,b);

 48:   /* Create PC context and set up data structures */
 49:   PCCreate(PETSC_COMM_WORLD,&pc);
 50:   PCSetType(pc,PCNONE);
 51:   PCSetFromOptions(pc);
 52:   PCSetOperators(pc,mat,mat);
 53:   PCSetUp(pc);

 55:   /* Create KSP context and set up data structures */
 56:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 57:   KSPSetType(ksp,KSPRICHARDSON);
 58:   KSPSetFromOptions(ksp);
 59:   PCSetOperators(pc,mat,mat);
 60:   KSPSetPC(ksp,pc);
 61:   KSPSetUp(ksp);

 63:   /* Solve the problem */
 64:   KSPGetType(ksp,&kspname);
 65:   PCGetType(pc,&pcname);
 66:   PetscPrintf(PETSC_COMM_SELF,"Running %s with %s preconditioning\n",kspname,pcname);
 67:   KSPSolve(ksp,b,u);
 68:   VecAXPY(u,-1.0,ustar);
 69:   VecNorm(u,NORM_2,&norm);
 70:   KSPGetIterationNumber(ksp,&its);
 71:   if (norm > tol) {
 72:     PetscPrintf(PETSC_COMM_SELF,"2 norm of error %g Number of iterations %D\n",(double)norm,its);
 73:   }

 75:   /* Free data structures */
 76:   KSPDestroy(&ksp);
 77:   VecDestroy(&u);
 78:   VecDestroy(&ustar);
 79:   VecDestroy(&b);
 80:   MatDestroy(&mat);
 81:   PCDestroy(&pc);

 83:   PetscFinalize();
 84:   return 0;
 85: }