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