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