xref: /petsc/src/ksp/pc/tutorials/ex1.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 static char help[] = "Test file for the PCFactorSetShiftType()\n";
2 /*
3  * Test file for the PCFactorSetShiftType() routine or -pc_factor_shift_type POSITIVE_DEFINITE option.
4  * The test matrix is the example from Kershaw's paper [J.Comp.Phys 1978]
5  * of a positive definite matrix for which ILU(0) will give a negative pivot.
6  * This means that the CG method will break down; the Manteuffel shift
7  * [Math. Comp. 1980] repairs this.
8  *
9  * Run the executable twice:
10  * 1/ without options: the iterative method diverges because of an
11  *    indefinite preconditioner
12  * 2/ with -pc_factor_shift_type POSITIVE_DEFINITE option (or comment in the PCFactorSetShiftType() line below):
13  *    the method will now successfully converge.
14  *
15  * Contributed by Victor Eijkhout 2003.
16  */
17 
18 #include <petscksp.h>
19 
20 int main(int argc,char **argv)
21 {
22   KSP                solver;
23   PC                 prec;
24   Mat                A,M;
25   Vec                X,B,D;
26   MPI_Comm           comm;
27   PetscScalar        v;
28   KSPConvergedReason reason;
29   PetscInt           i,j,its;
30   PetscErrorCode     ierr;
31 
32   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
33   comm = MPI_COMM_SELF;
34 
35   /*
36    * Construct the Kershaw matrix
37    * and a suitable rhs / initial guess
38    */
39   CHKERRQ(MatCreateSeqAIJ(comm,4,4,4,0,&A));
40   CHKERRQ(VecCreateSeq(comm,4,&B));
41   CHKERRQ(VecDuplicate(B,&X));
42   for (i=0; i<4; i++) {
43     v    = 3;
44     CHKERRQ(MatSetValues(A,1,&i,1,&i,&v,INSERT_VALUES));
45     v    = 1;
46     CHKERRQ(VecSetValues(B,1,&i,&v,INSERT_VALUES));
47     CHKERRQ(VecSetValues(X,1,&i,&v,INSERT_VALUES));
48   }
49 
50   i=0; v=0;
51   CHKERRQ(VecSetValues(X,1,&i,&v,INSERT_VALUES));
52 
53   for (i=0; i<3; i++) {
54     v    = -2; j=i+1;
55     CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
56     CHKERRQ(MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES));
57   }
58   i=0; j=3; v=2;
59 
60   CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
61   CHKERRQ(MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES));
62   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
63   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
64   CHKERRQ(VecAssemblyBegin(B));
65   CHKERRQ(VecAssemblyEnd(B));
66 
67   /*
68    * A Conjugate Gradient method
69    * with ILU(0) preconditioning
70    */
71   CHKERRQ(KSPCreate(comm,&solver));
72   CHKERRQ(KSPSetOperators(solver,A,A));
73 
74   CHKERRQ(KSPSetType(solver,KSPCG));
75   CHKERRQ(KSPSetInitialGuessNonzero(solver,PETSC_TRUE));
76 
77   /*
78    * ILU preconditioner;
79    * this will break down unless you add the Shift line,
80    * or use the -pc_factor_shift_positive_definite option */
81   CHKERRQ(KSPGetPC(solver,&prec));
82   CHKERRQ(PCSetType(prec,PCILU));
83   /* CHKERRQ(PCFactorSetShiftType(prec,MAT_SHIFT_POSITIVE_DEFINITE)); */
84 
85   CHKERRQ(KSPSetFromOptions(solver));
86   CHKERRQ(KSPSetUp(solver));
87 
88   /*
89    * Now that the factorisation is done, show the pivots;
90    * note that the last one is negative. This in itself is not an error,
91    * but it will make the iterative method diverge.
92    */
93   CHKERRQ(PCFactorGetMatrix(prec,&M));
94   CHKERRQ(VecDuplicate(B,&D));
95   CHKERRQ(MatGetDiagonal(M,D));
96 
97   /*
98    * Solve the system;
99    * without the shift this will diverge with
100    * an indefinite preconditioner
101    */
102   CHKERRQ(KSPSolve(solver,B,X));
103   CHKERRQ(KSPGetConvergedReason(solver,&reason));
104   if (reason==KSP_DIVERGED_INDEFINITE_PC) {
105     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nDivergence because of indefinite preconditioner;\n"));
106     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Run the executable again but with '-pc_factor_shift_type POSITIVE_DEFINITE' option.\n"));
107   } else if (reason<0) {
108     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nOther kind of divergence: this should not happen.\n"));
109   } else {
110     CHKERRQ(KSPGetIterationNumber(solver,&its));
111   }
112 
113   CHKERRQ(VecDestroy(&X));
114   CHKERRQ(VecDestroy(&B));
115   CHKERRQ(VecDestroy(&D));
116   CHKERRQ(MatDestroy(&A));
117   CHKERRQ(KSPDestroy(&solver));
118   ierr = PetscFinalize();
119   return ierr;
120 }
121 
122 /*TEST
123 
124    test:
125       args: -pc_factor_shift_type positive_definite
126 
127 TEST*/
128