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