Actual source code: ex1.c
petsc-3.6.1 2015-08-06
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: */
18: #include <petscksp.h>
22: int main(int argc,char **argv)
23: {
24: KSP solver;
25: PC prec;
26: Mat A,M;
27: Vec X,B,D;
28: MPI_Comm comm;
29: PetscScalar v;
30: KSPConvergedReason reason;
31: PetscInt i,j,its;
32: PetscErrorCode ierr;
34: PetscInitialize(&argc,&argv,0,help);
35: comm = MPI_COMM_SELF;
37: /*
38: * Construct the Kershaw matrix
39: * and a suitable rhs / initial guess
40: */
41: MatCreateSeqAIJ(comm,4,4,4,0,&A);
42: VecCreateSeq(comm,4,&B);
43: VecDuplicate(B,&X);
44: for (i=0; i<4; i++) {
45: v = 3;
46: MatSetValues(A,1,&i,1,&i,&v,INSERT_VALUES);
47: v = 1;
48: VecSetValues(B,1,&i,&v,INSERT_VALUES);
49: VecSetValues(X,1,&i,&v,INSERT_VALUES);
50: }
52: i=0; v=0;
53: VecSetValues(X,1,&i,&v,INSERT_VALUES);
55: for (i=0; i<3; i++) {
56: v = -2; j=i+1;
57: MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);
58: MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES);
59: }
60: i=0; j=3; v=2;
62: MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);
63: MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES);
64: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
65: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
66: VecAssemblyBegin(B);
67: VecAssemblyEnd(B);
69: /*
70: * A Conjugate Gradient method
71: * with ILU(0) preconditioning
72: */
73: KSPCreate(comm,&solver);
74: KSPSetOperators(solver,A,A);
76: KSPSetType(solver,KSPCG);
77: KSPSetInitialGuessNonzero(solver,PETSC_TRUE);
79: /*
80: * ILU preconditioner;
81: * this will break down unless you add the Shift line,
82: * or use the -pc_factor_shift_positive_definite option */
83: KSPGetPC(solver,&prec);
84: PCSetType(prec,PCILU);
85: /* PCFactorSetShiftType(prec,MAT_SHIFT_POSITIVE_DEFINITE); */
87: KSPSetFromOptions(solver);
88: KSPSetUp(solver);
90: /*
91: * Now that the factorisation is done, show the pivots;
92: * note that the last one is negative. This in itself is not an error,
93: * but it will make the iterative method diverge.
94: */
95: PCFactorGetMatrix(prec,&M);
96: VecDuplicate(B,&D);
97: MatGetDiagonal(M,D);
99: /*
100: * Solve the system;
101: * without the shift this will diverge with
102: * an indefinite preconditioner
103: */
104: KSPSolve(solver,B,X);
105: KSPGetConvergedReason(solver,&reason);
106: if (reason==KSP_DIVERGED_INDEFINITE_PC) {
107: PetscPrintf(PETSC_COMM_WORLD,"\nDivergence because of indefinite preconditioner;\n");
108: PetscPrintf(PETSC_COMM_WORLD,"Run the executable again but with '-pc_factor_shift_type POSITIVE_DEFINITE' option.\n");
109: } else if (reason<0) {
110: PetscPrintf(PETSC_COMM_WORLD,"\nOther kind of divergence: this should not happen.\n");
111: } else {
112: KSPGetIterationNumber(solver,&its);
113: }
115: VecDestroy(&X);
116: VecDestroy(&B);
117: VecDestroy(&D);
118: MatDestroy(&A);
119: KSPDestroy(&solver);
120: PetscFinalize();
121: return 0;
122: }