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 18c4762a1bSJed Brown int main(int argc,char **argv) 19c4762a1bSJed Brown { 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; 30*327415f7SBarry 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 49c4762a1bSJed Brown i =0; v=0; 509566063dSJacob Faibussowitsch PetscCall(VecSetValues(X,1,&i,&v,INSERT_VALUES)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown for (i=0; i<3; i++) { 53c4762a1bSJed Brown v = -2; j=i+1; 549566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 559566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES)); 56c4762a1bSJed Brown } 57c4762a1bSJed Brown i=0; j=3; v=2; 58c4762a1bSJed Brown 599566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 609566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES)); 619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 639566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(B)); 649566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(B)); 659566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nThe Kershaw matrix:\n\n")); 669566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* 69c4762a1bSJed Brown * A Conjugate Gradient method 70c4762a1bSJed Brown * with ILU(0) preconditioning 71c4762a1bSJed Brown */ 729566063dSJacob Faibussowitsch PetscCall(KSPCreate(comm,&ksp)); 739566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp,A,A)); 74c4762a1bSJed Brown 759566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp,KSPCG)); 769566063dSJacob Faibussowitsch PetscCall(KSPSetInitialGuessNonzero(ksp,PETSC_TRUE)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* 79c4762a1bSJed Brown * ILU preconditioner; 80c4762a1bSJed Brown * The iterative method will break down unless you comment in the SetShift 81c4762a1bSJed Brown * line below, or use the -pc_factor_shift_positive_definite option. 82c4762a1bSJed Brown * Run the code twice: once as given to see the negative pivot and the 83c4762a1bSJed Brown * divergence behaviour, then comment in the Shift line, or add the 84c4762a1bSJed Brown * command line option, and see that the pivots are all positive and 85c4762a1bSJed Brown * the method converges. 86c4762a1bSJed Brown */ 879566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 889566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCICC)); 899566063dSJacob Faibussowitsch /* PetscCall(PCFactorSetShiftType(prec,MAT_SHIFT_POSITIVE_DEFINITE)); */ 90c4762a1bSJed Brown 919566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 929566063dSJacob Faibussowitsch PetscCall(KSPSetUp(ksp)); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* 95c4762a1bSJed Brown * Now that the factorisation is done, show the pivots; 96c4762a1bSJed Brown * note that the last one is negative. This in itself is not an error, 97c4762a1bSJed Brown * but it will make the iterative method diverge. 98c4762a1bSJed Brown */ 999566063dSJacob Faibussowitsch PetscCall(PCFactorGetMatrix(pc,&M)); 1009566063dSJacob Faibussowitsch PetscCall(VecDuplicate(B,&D)); 1019566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(M,D)); 1029566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nPivots:\n\n")); 1039566063dSJacob Faibussowitsch PetscCall(VecView(D,0)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* 106c4762a1bSJed Brown * Solve the system; 107c4762a1bSJed Brown * without the shift this will diverge with 108c4762a1bSJed Brown * an indefinite preconditioner 109c4762a1bSJed Brown */ 1109566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp,B,X)); 1119566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(ksp,&reason)); 112c4762a1bSJed Brown if (reason==KSP_DIVERGED_INDEFINITE_PC) { 1139566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nDivergence because of indefinite preconditioner;\n")); 1149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Run the executable again but with -pc_factor_shift_positive_definite option.\n")); 115c4762a1bSJed Brown } else if (reason<0) { 1169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nOther kind of divergence: this should not happen.\n")); 117c4762a1bSJed Brown } else { 1189566063dSJacob Faibussowitsch PetscCall(KSPGetIterationNumber(ksp,&its)); 1199566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nConvergence in %d iterations.\n",(int)its)); 120c4762a1bSJed Brown } 1219566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n")); 122c4762a1bSJed Brown 1239566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 1249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&B)); 1269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 1279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&D)); 1289566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 129b122ec5aSJacob Faibussowitsch return 0; 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown /*TEST 133c4762a1bSJed Brown 134c4762a1bSJed Brown test: 135560a203cSprj- filter: sed -e "s/in 5 iterations/in 4 iterations/g" 136560a203cSprj- 137c4762a1bSJed Brown TEST*/ 138