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