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
main(int argc,char ** argv)20d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
21d71ae5a4SJacob Faibussowitsch {
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
31327415f7SBarry Smith PetscFunctionBeginUser;
329566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help));
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 */
399566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(comm, 4, 4, 4, 0, &A));
409566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(comm, 4, &B));
419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(B, &X));
42c4762a1bSJed Brown for (i = 0; i < 4; i++) {
43c4762a1bSJed Brown v = 3;
449566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &i, &v, INSERT_VALUES));
45c4762a1bSJed Brown v = 1;
469566063dSJacob Faibussowitsch PetscCall(VecSetValues(B, 1, &i, &v, INSERT_VALUES));
479566063dSJacob Faibussowitsch PetscCall(VecSetValues(X, 1, &i, &v, INSERT_VALUES));
48c4762a1bSJed Brown }
49c4762a1bSJed Brown
509371c9d4SSatish Balay i = 0;
519371c9d4SSatish Balay v = 0;
529566063dSJacob Faibussowitsch PetscCall(VecSetValues(X, 1, &i, &v, INSERT_VALUES));
53c4762a1bSJed Brown
54c4762a1bSJed Brown for (i = 0; i < 3; i++) {
559371c9d4SSatish Balay v = -2;
569371c9d4SSatish Balay j = i + 1;
579566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
589566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &j, 1, &i, &v, INSERT_VALUES));
59c4762a1bSJed Brown }
609371c9d4SSatish Balay i = 0;
619371c9d4SSatish Balay j = 3;
629371c9d4SSatish Balay v = 2;
63c4762a1bSJed Brown
649566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES));
659566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &j, 1, &i, &v, INSERT_VALUES));
669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
689566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(B));
699566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(B));
70c4762a1bSJed Brown
71c4762a1bSJed Brown /*
72c4762a1bSJed Brown * A Conjugate Gradient method
73c4762a1bSJed Brown * with ILU(0) preconditioning
74c4762a1bSJed Brown */
759566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &solver));
769566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(solver, A, A));
77c4762a1bSJed Brown
789566063dSJacob Faibussowitsch PetscCall(KSPSetType(solver, KSPCG));
799566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(solver, PETSC_TRUE));
80c4762a1bSJed Brown
81c4762a1bSJed Brown /*
82c4762a1bSJed Brown * ILU preconditioner;
83c4762a1bSJed Brown * this will break down unless you add the Shift line,
84c4762a1bSJed Brown * or use the -pc_factor_shift_positive_definite option */
859566063dSJacob Faibussowitsch PetscCall(KSPGetPC(solver, &prec));
869566063dSJacob Faibussowitsch PetscCall(PCSetType(prec, PCILU));
879566063dSJacob Faibussowitsch /* PetscCall(PCFactorSetShiftType(prec,MAT_SHIFT_POSITIVE_DEFINITE)); */
88c4762a1bSJed Brown
899566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(solver));
909566063dSJacob Faibussowitsch PetscCall(KSPSetUp(solver));
91c4762a1bSJed Brown
92c4762a1bSJed Brown /*
93c4762a1bSJed Brown * Now that the factorisation is done, show the pivots;
94c4762a1bSJed Brown * note that the last one is negative. This in itself is not an error,
95c4762a1bSJed Brown * but it will make the iterative method diverge.
96c4762a1bSJed Brown */
979566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(prec, &M));
989566063dSJacob Faibussowitsch PetscCall(VecDuplicate(B, &D));
999566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(M, D));
100c4762a1bSJed Brown
101c4762a1bSJed Brown /*
102c4762a1bSJed Brown * Solve the system;
103c4762a1bSJed Brown * without the shift this will diverge with
104c4762a1bSJed Brown * an indefinite preconditioner
105c4762a1bSJed Brown */
1069566063dSJacob Faibussowitsch PetscCall(KSPSolve(solver, B, X));
1079566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(solver, &reason));
108c4762a1bSJed Brown if (reason == KSP_DIVERGED_INDEFINITE_PC) {
1099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nDivergence because of indefinite preconditioner;\n"));
1109566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Run the executable again but with '-pc_factor_shift_type POSITIVE_DEFINITE' option.\n"));
111c4762a1bSJed Brown } else if (reason < 0) {
1129566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nOther kind of divergence: this should not happen.\n"));
113c4762a1bSJed Brown } else {
1149566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(solver, &its));
115c4762a1bSJed Brown }
116c4762a1bSJed Brown
1179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
1189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&B));
1199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&D));
1209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1219566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&solver));
1229566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
123b122ec5aSJacob Faibussowitsch return 0;
124c4762a1bSJed Brown }
125c4762a1bSJed Brown
126c4762a1bSJed Brown /*TEST
127c4762a1bSJed Brown
128c4762a1bSJed Brown test:
129c4762a1bSJed Brown args: -pc_factor_shift_type positive_definite
130*3886731fSPierre Jolivet output_file: output/empty.out
131c4762a1bSJed Brown
132c4762a1bSJed Brown TEST*/
133