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 18*9371c9d4SSatish Balay int main(int argc, char **argv) { 19c4762a1bSJed Brown KSP ksp; 20c4762a1bSJed Brown PC pc; 21c4762a1bSJed Brown Mat A, M; 22c4762a1bSJed Brown Vec X, B, D; 23c4762a1bSJed Brown MPI_Comm comm; 24c4762a1bSJed Brown PetscScalar v; 25c4762a1bSJed Brown KSPConvergedReason reason; 26c4762a1bSJed Brown PetscInt i, j, its; 27c4762a1bSJed Brown 28c4762a1bSJed Brown PetscFunctionBegin; 29327415f7SBarry Smith PetscFunctionBeginUser; 309566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 31c4762a1bSJed Brown comm = MPI_COMM_SELF; 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* 34c4762a1bSJed Brown * Construct the Kershaw matrix 35c4762a1bSJed Brown * and a suitable rhs / initial guess 36c4762a1bSJed Brown */ 379566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(comm, 4, 4, 4, 0, &A)); 389566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(comm, 4, &B)); 399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(B, &X)); 40c4762a1bSJed Brown for (i = 0; i < 4; i++) { 41c4762a1bSJed Brown v = 3; 429566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &i, &v, INSERT_VALUES)); 43c4762a1bSJed Brown v = 1; 449566063dSJacob Faibussowitsch PetscCall(VecSetValues(B, 1, &i, &v, INSERT_VALUES)); 459566063dSJacob Faibussowitsch PetscCall(VecSetValues(X, 1, &i, &v, INSERT_VALUES)); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48*9371c9d4SSatish Balay i = 0; 49*9371c9d4SSatish Balay v = 0; 509566063dSJacob Faibussowitsch PetscCall(VecSetValues(X, 1, &i, &v, INSERT_VALUES)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown for (i = 0; i < 3; i++) { 53*9371c9d4SSatish Balay v = -2; 54*9371c9d4SSatish Balay j = i + 1; 559566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 569566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &j, 1, &i, &v, INSERT_VALUES)); 57c4762a1bSJed Brown } 58*9371c9d4SSatish Balay i = 0; 59*9371c9d4SSatish Balay j = 3; 60*9371c9d4SSatish Balay v = 2; 61c4762a1bSJed Brown 629566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &i, 1, &j, &v, INSERT_VALUES)); 639566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &j, 1, &i, &v, INSERT_VALUES)); 649566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 659566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 669566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(B)); 679566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(B)); 689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nThe Kershaw matrix:\n\n")); 699566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown /* 72c4762a1bSJed Brown * A Conjugate Gradient method 73c4762a1bSJed Brown * with ILU(0) preconditioning 74c4762a1bSJed Brown */ 759566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm, &ksp)); 769566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp, A, A)); 77c4762a1bSJed Brown 789566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp, KSPCG)); 799566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(ksp, PETSC_TRUE)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* 82c4762a1bSJed Brown * ILU preconditioner; 83c4762a1bSJed Brown * The iterative method will break down unless you comment in the SetShift 84c4762a1bSJed Brown * line below, or use the -pc_factor_shift_positive_definite option. 85c4762a1bSJed Brown * Run the code twice: once as given to see the negative pivot and the 86c4762a1bSJed Brown * divergence behaviour, then comment in the Shift line, or add the 87c4762a1bSJed Brown * command line option, and see that the pivots are all positive and 88c4762a1bSJed Brown * the method converges. 89c4762a1bSJed Brown */ 909566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 919566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCICC)); 929566063dSJacob Faibussowitsch /* PetscCall(PCFactorSetShiftType(prec,MAT_SHIFT_POSITIVE_DEFINITE)); */ 93c4762a1bSJed Brown 949566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 959566063dSJacob Faibussowitsch PetscCall(KSPSetUp(ksp)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* 98c4762a1bSJed Brown * Now that the factorisation is done, show the pivots; 99c4762a1bSJed Brown * note that the last one is negative. This in itself is not an error, 100c4762a1bSJed Brown * but it will make the iterative method diverge. 101c4762a1bSJed Brown */ 1029566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc, &M)); 1039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(B, &D)); 1049566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(M, D)); 1059566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nPivots:\n\n")); 1069566063dSJacob Faibussowitsch PetscCall(VecView(D, 0)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* 109c4762a1bSJed Brown * Solve the system; 110c4762a1bSJed Brown * without the shift this will diverge with 111c4762a1bSJed Brown * an indefinite preconditioner 112c4762a1bSJed Brown */ 1139566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp, B, X)); 1149566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp, &reason)); 115c4762a1bSJed Brown if (reason == KSP_DIVERGED_INDEFINITE_PC) { 1169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nDivergence because of indefinite preconditioner;\n")); 1179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Run the executable again but with -pc_factor_shift_positive_definite option.\n")); 118c4762a1bSJed Brown } else if (reason < 0) { 1199566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nOther kind of divergence: this should not happen.\n")); 120c4762a1bSJed Brown } else { 1219566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(ksp, &its)); 1229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nConvergence in %d iterations.\n", (int)its)); 123c4762a1bSJed Brown } 1249566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n")); 125c4762a1bSJed Brown 1269566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 1279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&B)); 1299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 1309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&D)); 1319566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 132b122ec5aSJacob Faibussowitsch return 0; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /*TEST 136c4762a1bSJed Brown 137c4762a1bSJed Brown test: 138560a203cSprj- filter: sed -e "s/in 5 iterations/in 4 iterations/g" 139560a203cSprj- 140c4762a1bSJed Brown TEST*/ 141