xref: /petsc/src/ksp/pc/tutorials/ex1.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode     ierr;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
33c4762a1bSJed Brown   comm = MPI_COMM_SELF;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   /*
36c4762a1bSJed Brown    * Construct the Kershaw matrix
37c4762a1bSJed Brown    * and a suitable rhs / initial guess
38c4762a1bSJed Brown    */
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(comm,4,4,4,0,&A));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(comm,4,&B));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(B,&X));
42c4762a1bSJed Brown   for (i=0; i<4; i++) {
43c4762a1bSJed Brown     v    = 3;
44*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&i,1,&i,&v,INSERT_VALUES));
45c4762a1bSJed Brown     v    = 1;
46*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(B,1,&i,&v,INSERT_VALUES));
47*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(X,1,&i,&v,INSERT_VALUES));
48c4762a1bSJed Brown   }
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   i=0; v=0;
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetValues(X,1,&i,&v,INSERT_VALUES));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   for (i=0; i<3; i++) {
54c4762a1bSJed Brown     v    = -2; j=i+1;
55*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
56*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES));
57c4762a1bSJed Brown   }
58c4762a1bSJed Brown   i=0; j=3; v=2;
59c4762a1bSJed Brown 
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES));
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(B));
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(B));
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   /*
68c4762a1bSJed Brown    * A Conjugate Gradient method
69c4762a1bSJed Brown    * with ILU(0) preconditioning
70c4762a1bSJed Brown    */
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPCreate(comm,&solver));
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetOperators(solver,A,A));
73c4762a1bSJed Brown 
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetType(solver,KSPCG));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetInitialGuessNonzero(solver,PETSC_TRUE));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /*
78c4762a1bSJed Brown    * ILU preconditioner;
79c4762a1bSJed Brown    * this will break down unless you add the Shift line,
80c4762a1bSJed Brown    * or use the -pc_factor_shift_positive_definite option */
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(solver,&prec));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCSetType(prec,PCILU));
83*5f80ce2aSJacob Faibussowitsch   /* CHKERRQ(PCFactorSetShiftType(prec,MAT_SHIFT_POSITIVE_DEFINITE)); */
84c4762a1bSJed Brown 
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetFromOptions(solver));
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSetUp(solver));
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   /*
89c4762a1bSJed Brown    * Now that the factorisation is done, show the pivots;
90c4762a1bSJed Brown    * note that the last one is negative. This in itself is not an error,
91c4762a1bSJed Brown    * but it will make the iterative method diverge.
92c4762a1bSJed Brown    */
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PCFactorGetMatrix(prec,&M));
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(B,&D));
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetDiagonal(M,D));
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /*
98c4762a1bSJed Brown    * Solve the system;
99c4762a1bSJed Brown    * without the shift this will diverge with
100c4762a1bSJed Brown    * an indefinite preconditioner
101c4762a1bSJed Brown    */
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPSolve(solver,B,X));
103*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetConvergedReason(solver,&reason));
104c4762a1bSJed Brown   if (reason==KSP_DIVERGED_INDEFINITE_PC) {
105*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nDivergence because of indefinite preconditioner;\n"));
106*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Run the executable again but with '-pc_factor_shift_type POSITIVE_DEFINITE' option.\n"));
107c4762a1bSJed Brown   } else if (reason<0) {
108*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nOther kind of divergence: this should not happen.\n"));
109c4762a1bSJed Brown   } else {
110*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetIterationNumber(solver,&its));
111c4762a1bSJed Brown   }
112c4762a1bSJed Brown 
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&B));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&D));
116*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPDestroy(&solver));
118c4762a1bSJed Brown   ierr = PetscFinalize();
119c4762a1bSJed Brown   return ierr;
120c4762a1bSJed Brown }
121c4762a1bSJed Brown 
122c4762a1bSJed Brown /*TEST
123c4762a1bSJed Brown 
124c4762a1bSJed Brown    test:
125c4762a1bSJed Brown       args: -pc_factor_shift_type positive_definite
126c4762a1bSJed Brown 
127c4762a1bSJed Brown TEST*/
128