xref: /petsc/src/ksp/pc/tutorials/ex1.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
20*9371c9d4SSatish Balay int main(int argc, char **argv) {
21c4762a1bSJed Brown   KSP                solver;
22c4762a1bSJed Brown   PC                 prec;
23c4762a1bSJed Brown   Mat                A, M;
24c4762a1bSJed Brown   Vec                X, B, D;
25c4762a1bSJed Brown   MPI_Comm           comm;
26c4762a1bSJed Brown   PetscScalar        v;
27c4762a1bSJed Brown   KSPConvergedReason reason;
28c4762a1bSJed Brown   PetscInt           i, j, its;
29c4762a1bSJed Brown 
30327415f7SBarry Smith   PetscFunctionBeginUser;
319566063dSJacob 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    */
389566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(comm, 4, 4, 4, 0, &A));
399566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(comm, 4, &B));
409566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(B, &X));
41c4762a1bSJed Brown   for (i = 0; i < 4; i++) {
42c4762a1bSJed Brown     v = 3;
439566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 1, &i, &v, INSERT_VALUES));
44c4762a1bSJed Brown     v = 1;
459566063dSJacob Faibussowitsch     PetscCall(VecSetValues(B, 1, &i, &v, INSERT_VALUES));
469566063dSJacob Faibussowitsch     PetscCall(VecSetValues(X, 1, &i, &v, INSERT_VALUES));
47c4762a1bSJed Brown   }
48c4762a1bSJed Brown 
49*9371c9d4SSatish Balay   i = 0;
50*9371c9d4SSatish Balay   v = 0;
519566063dSJacob Faibussowitsch   PetscCall(VecSetValues(X, 1, &i, &v, INSERT_VALUES));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
54*9371c9d4SSatish Balay     v = -2;
55*9371c9d4SSatish Balay     j = i + 1;
569566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
579566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &j, 1, &i, &v, INSERT_VALUES));
58c4762a1bSJed Brown   }
59*9371c9d4SSatish Balay   i = 0;
60*9371c9d4SSatish Balay   j = 3;
61*9371c9d4SSatish Balay   v = 2;
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
649566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &j, 1, &i, &v, INSERT_VALUES));
659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
679566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(B));
689566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(B));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /*
71c4762a1bSJed Brown    * A Conjugate Gradient method
72c4762a1bSJed Brown    * with ILU(0) preconditioning
73c4762a1bSJed Brown    */
749566063dSJacob Faibussowitsch   PetscCall(KSPCreate(comm, &solver));
759566063dSJacob Faibussowitsch   PetscCall(KSPSetOperators(solver, A, A));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(KSPSetType(solver, KSPCG));
789566063dSJacob Faibussowitsch   PetscCall(KSPSetInitialGuessNonzero(solver, PETSC_TRUE));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /*
81c4762a1bSJed Brown    * ILU preconditioner;
82c4762a1bSJed Brown    * this will break down unless you add the Shift line,
83c4762a1bSJed Brown    * or use the -pc_factor_shift_positive_definite option */
849566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(solver, &prec));
859566063dSJacob Faibussowitsch   PetscCall(PCSetType(prec, PCILU));
869566063dSJacob Faibussowitsch   /* PetscCall(PCFactorSetShiftType(prec,MAT_SHIFT_POSITIVE_DEFINITE)); */
87c4762a1bSJed Brown 
889566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(solver));
899566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(solver));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /*
92c4762a1bSJed Brown    * Now that the factorisation is done, show the pivots;
93c4762a1bSJed Brown    * note that the last one is negative. This in itself is not an error,
94c4762a1bSJed Brown    * but it will make the iterative method diverge.
95c4762a1bSJed Brown    */
969566063dSJacob Faibussowitsch   PetscCall(PCFactorGetMatrix(prec, &M));
979566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(B, &D));
989566063dSJacob Faibussowitsch   PetscCall(MatGetDiagonal(M, D));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /*
101c4762a1bSJed Brown    * Solve the system;
102c4762a1bSJed Brown    * without the shift this will diverge with
103c4762a1bSJed Brown    * an indefinite preconditioner
104c4762a1bSJed Brown    */
1059566063dSJacob Faibussowitsch   PetscCall(KSPSolve(solver, B, X));
1069566063dSJacob Faibussowitsch   PetscCall(KSPGetConvergedReason(solver, &reason));
107c4762a1bSJed Brown   if (reason == KSP_DIVERGED_INDEFINITE_PC) {
1089566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nDivergence because of indefinite preconditioner;\n"));
1099566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Run the executable again but with '-pc_factor_shift_type POSITIVE_DEFINITE' option.\n"));
110c4762a1bSJed Brown   } else if (reason < 0) {
1119566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nOther kind of divergence: this should not happen.\n"));
112c4762a1bSJed Brown   } else {
1139566063dSJacob Faibussowitsch     PetscCall(KSPGetIterationNumber(solver, &its));
114c4762a1bSJed Brown   }
115c4762a1bSJed Brown 
1169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&B));
1189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&D));
1199566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1209566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&solver));
1219566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
122b122ec5aSJacob Faibussowitsch   return 0;
123c4762a1bSJed Brown }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown /*TEST
126c4762a1bSJed Brown 
127c4762a1bSJed Brown    test:
128c4762a1bSJed Brown       args: -pc_factor_shift_type positive_definite
129c4762a1bSJed Brown 
130c4762a1bSJed Brown TEST*/
131