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_positive_definite option (or comment in the PCFactorSetShiftType() line below):
13c4762a1bSJed Brown * the method will now successfully converge.
14c4762a1bSJed Brown */
15c4762a1bSJed Brown
16c4762a1bSJed Brown #include <petscksp.h>
17c4762a1bSJed Brown
main(int argc,char ** argv)18d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
19d71ae5a4SJacob Faibussowitsch {
20c4762a1bSJed Brown KSP ksp;
21c4762a1bSJed Brown PC pc;
22c4762a1bSJed Brown Mat A, M;
23c4762a1bSJed Brown Vec X, B, D;
24c4762a1bSJed Brown MPI_Comm comm;
25c4762a1bSJed Brown PetscScalar v;
26c4762a1bSJed Brown KSPConvergedReason reason;
27c4762a1bSJed Brown PetscInt i, j, its;
28c4762a1bSJed Brown
29c4762a1bSJed Brown PetscFunctionBegin;
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
499371c9d4SSatish Balay i = 0;
509371c9d4SSatish Balay v = 0;
519566063dSJacob Faibussowitsch PetscCall(VecSetValues(X, 1, &i, &v, INSERT_VALUES));
52c4762a1bSJed Brown
53c4762a1bSJed Brown for (i = 0; i < 3; i++) {
549371c9d4SSatish Balay v = -2;
559371c9d4SSatish 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 }
599371c9d4SSatish Balay i = 0;
609371c9d4SSatish Balay j = 3;
619371c9d4SSatish 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));
699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nThe Kershaw matrix:\n\n"));
709566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
71c4762a1bSJed Brown
72c4762a1bSJed Brown /*
73c4762a1bSJed Brown * A Conjugate Gradient method
74c4762a1bSJed Brown * with ILU(0) preconditioning
75c4762a1bSJed Brown */
769566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp));
779566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, A, A));
78c4762a1bSJed Brown
799566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp, KSPCG));
809566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_TRUE));
81c4762a1bSJed Brown
82c4762a1bSJed Brown /*
83c4762a1bSJed Brown * ILU preconditioner;
84c4762a1bSJed Brown * The iterative method will break down unless you comment in the SetShift
85c4762a1bSJed Brown * line below, or use the -pc_factor_shift_positive_definite option.
86c4762a1bSJed Brown * Run the code twice: once as given to see the negative pivot and the
87c4762a1bSJed Brown * divergence behaviour, then comment in the Shift line, or add the
88c4762a1bSJed Brown * command line option, and see that the pivots are all positive and
89c4762a1bSJed Brown * the method converges.
90c4762a1bSJed Brown */
919566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc));
929566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCICC));
939566063dSJacob Faibussowitsch /* PetscCall(PCFactorSetShiftType(prec,MAT_SHIFT_POSITIVE_DEFINITE)); */
94c4762a1bSJed Brown
959566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp));
969566063dSJacob Faibussowitsch PetscCall(KSPSetUp(ksp));
97c4762a1bSJed Brown
98c4762a1bSJed Brown /*
99c4762a1bSJed Brown * Now that the factorisation is done, show the pivots;
100c4762a1bSJed Brown * note that the last one is negative. This in itself is not an error,
101c4762a1bSJed Brown * but it will make the iterative method diverge.
102c4762a1bSJed Brown */
1039566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc, &M));
1049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(B, &D));
1059566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(M, D));
1069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nPivots:\n\n"));
1079566063dSJacob Faibussowitsch PetscCall(VecView(D, 0));
108c4762a1bSJed Brown
109c4762a1bSJed Brown /*
110c4762a1bSJed Brown * Solve the system;
111c4762a1bSJed Brown * without the shift this will diverge with
112c4762a1bSJed Brown * an indefinite preconditioner
113c4762a1bSJed Brown */
1149566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, B, X));
1159566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &reason));
116c4762a1bSJed Brown if (reason == KSP_DIVERGED_INDEFINITE_PC) {
1179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nDivergence because of indefinite preconditioner;\n"));
1189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Run the executable again but with -pc_factor_shift_positive_definite option.\n"));
119c4762a1bSJed Brown } else if (reason < 0) {
1209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nOther kind of divergence: this should not happen.\n"));
121c4762a1bSJed Brown } else {
1229566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(ksp, &its));
123*300f1712SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nConvergence in %" PetscInt_FMT " iterations.\n", its));
124c4762a1bSJed Brown }
1259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
126c4762a1bSJed Brown
1279566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp));
1289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
1299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&B));
1309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
1319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&D));
1329566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
133b122ec5aSJacob Faibussowitsch return 0;
134c4762a1bSJed Brown }
135c4762a1bSJed Brown
136c4762a1bSJed Brown /*TEST
137c4762a1bSJed Brown
138c4762a1bSJed Brown test:
139560a203cSprj- filter: sed -e "s/in 5 iterations/in 4 iterations/g"
140560a203cSprj-
141c4762a1bSJed Brown TEST*/
142