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*b122ec5aSJacob Faibussowitsch CHKERRQ(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 */ 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqAIJ(comm,4,4,4,0,&A)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(comm,4,&B)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(B,&X)); 40c4762a1bSJed Brown for (i=0; i<4; i++) { 41c4762a1bSJed Brown v = 3; 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&i,&v,INSERT_VALUES)); 43c4762a1bSJed Brown v = 1; 445f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(B,1,&i,&v,INSERT_VALUES)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(X,1,&i,&v,INSERT_VALUES)); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48c4762a1bSJed Brown i =0; v=0; 495f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(X,1,&i,&v,INSERT_VALUES)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown for (i=0; i<3; i++) { 52c4762a1bSJed Brown v = -2; j=i+1; 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES)); 55c4762a1bSJed Brown } 56c4762a1bSJed Brown i=0; j=3; v=2; 57c4762a1bSJed Brown 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(B)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(B)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nThe Kershaw matrix:\n\n")); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* 68c4762a1bSJed Brown * A Conjugate Gradient method 69c4762a1bSJed Brown * with ILU(0) preconditioning 70c4762a1bSJed Brown */ 715f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCreate(comm,&ksp)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(ksp,A,A)); 73c4762a1bSJed Brown 745f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetType(ksp,KSPCG)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetInitialGuessNonzero(ksp,PETSC_TRUE)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* 78c4762a1bSJed Brown * ILU preconditioner; 79c4762a1bSJed Brown * The iterative method will break down unless you comment in the SetShift 80c4762a1bSJed Brown * line below, or use the -pc_factor_shift_positive_definite option. 81c4762a1bSJed Brown * Run the code twice: once as given to see the negative pivot and the 82c4762a1bSJed Brown * divergence behaviour, then comment in the Shift line, or add the 83c4762a1bSJed Brown * command line option, and see that the pivots are all positive and 84c4762a1bSJed Brown * the method converges. 85c4762a1bSJed Brown */ 865f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(ksp,&pc)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pc,PCICC)); 885f80ce2aSJacob Faibussowitsch /* CHKERRQ(PCFactorSetShiftType(prec,MAT_SHIFT_POSITIVE_DEFINITE)); */ 89c4762a1bSJed Brown 905f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetFromOptions(ksp)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetUp(ksp)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* 94c4762a1bSJed Brown * Now that the factorisation is done, show the pivots; 95c4762a1bSJed Brown * note that the last one is negative. This in itself is not an error, 96c4762a1bSJed Brown * but it will make the iterative method diverge. 97c4762a1bSJed Brown */ 985f80ce2aSJacob Faibussowitsch CHKERRQ(PCFactorGetMatrix(pc,&M)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(B,&D)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(M,D)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nPivots:\n\n")); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(D,0)); 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* 105c4762a1bSJed Brown * Solve the system; 106c4762a1bSJed Brown * without the shift this will diverge with 107c4762a1bSJed Brown * an indefinite preconditioner 108c4762a1bSJed Brown */ 1095f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolve(ksp,B,X)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetConvergedReason(ksp,&reason)); 111c4762a1bSJed Brown if (reason==KSP_DIVERGED_INDEFINITE_PC) { 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nDivergence because of indefinite preconditioner;\n")); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Run the executable again but with -pc_factor_shift_positive_definite option.\n")); 114c4762a1bSJed Brown } else if (reason<0) { 1155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nOther kind of divergence: this should not happen.\n")); 116c4762a1bSJed Brown } else { 1175f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetIterationNumber(ksp,&its)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nConvergence in %d iterations.\n",(int)its)); 119c4762a1bSJed Brown } 1205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n")); 121c4762a1bSJed Brown 1225f80ce2aSJacob Faibussowitsch CHKERRQ(KSPDestroy(&ksp)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&B)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&D)); 127*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 128*b122ec5aSJacob Faibussowitsch return 0; 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 131c4762a1bSJed Brown /*TEST 132c4762a1bSJed Brown 133c4762a1bSJed Brown test: 134560a203cSprj- filter: sed -e "s/in 5 iterations/in 4 iterations/g" 135560a203cSprj- 136c4762a1bSJed Brown TEST*/ 137