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