1*c4762a1bSJed Brown static char help[] = "Test file for the PCFactorSetShiftType()\n"; 2*c4762a1bSJed Brown /* 3*c4762a1bSJed Brown * Test file for the PCFactorSetShiftType() routine or -pc_factor_shift_type POSITIVE_DEFINITE option. 4*c4762a1bSJed Brown * The test matrix is the example from Kershaw's paper [J.Comp.Phys 1978] 5*c4762a1bSJed Brown * of a positive definite matrix for which ILU(0) will give a negative pivot. 6*c4762a1bSJed Brown * This means that the CG method will break down; the Manteuffel shift 7*c4762a1bSJed Brown * [Math. Comp. 1980] repairs this. 8*c4762a1bSJed Brown * 9*c4762a1bSJed Brown * Run the executable twice: 10*c4762a1bSJed Brown * 1/ without options: the iterative method diverges because of an 11*c4762a1bSJed Brown * indefinite preconditioner 12*c4762a1bSJed Brown * 2/ with -pc_factor_shift_positive_definite option (or comment in the PCFactorSetShiftType() line below): 13*c4762a1bSJed Brown * the method will now successfully converge. 14*c4762a1bSJed Brown */ 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown #include <petscksp.h> 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown int main(int argc,char **argv) 19*c4762a1bSJed Brown { 20*c4762a1bSJed Brown KSP ksp; 21*c4762a1bSJed Brown PC pc; 22*c4762a1bSJed Brown Mat A,M; 23*c4762a1bSJed Brown Vec X,B,D; 24*c4762a1bSJed Brown MPI_Comm comm; 25*c4762a1bSJed Brown PetscScalar v; 26*c4762a1bSJed Brown KSPConvergedReason reason; 27*c4762a1bSJed Brown PetscInt i,j,its; 28*c4762a1bSJed Brown PetscErrorCode ierr; 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown PetscFunctionBegin; 31*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 32*c4762a1bSJed Brown comm = MPI_COMM_SELF; 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown /* 35*c4762a1bSJed Brown * Construct the Kershaw matrix 36*c4762a1bSJed Brown * and a suitable rhs / initial guess 37*c4762a1bSJed Brown */ 38*c4762a1bSJed Brown ierr = MatCreateSeqAIJ(comm,4,4,4,0,&A);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = VecCreateSeq(comm,4,&B);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = VecDuplicate(B,&X);CHKERRQ(ierr); 41*c4762a1bSJed Brown for (i=0; i<4; i++) { 42*c4762a1bSJed Brown v = 3; 43*c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 44*c4762a1bSJed Brown v = 1; 45*c4762a1bSJed Brown ierr = VecSetValues(B,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 46*c4762a1bSJed Brown ierr = VecSetValues(X,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 47*c4762a1bSJed Brown } 48*c4762a1bSJed Brown 49*c4762a1bSJed Brown i =0; v=0; 50*c4762a1bSJed Brown ierr = VecSetValues(X,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown for (i=0; i<3; i++) { 53*c4762a1bSJed Brown v = -2; j=i+1; 54*c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 56*c4762a1bSJed Brown } 57*c4762a1bSJed Brown i=0; j=3; v=2; 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = VecAssemblyBegin(B);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = VecAssemblyEnd(B);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nThe Kershaw matrix:\n\n");CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 67*c4762a1bSJed Brown 68*c4762a1bSJed Brown /* 69*c4762a1bSJed Brown * A Conjugate Gradient method 70*c4762a1bSJed Brown * with ILU(0) preconditioning 71*c4762a1bSJed Brown */ 72*c4762a1bSJed Brown ierr = KSPCreate(comm,&ksp);CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = KSPSetOperators(ksp,A,A);CHKERRQ(ierr); 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown ierr = KSPSetType(ksp,KSPCG);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr); 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown /* 79*c4762a1bSJed Brown * ILU preconditioner; 80*c4762a1bSJed Brown * The iterative method will break down unless you comment in the SetShift 81*c4762a1bSJed Brown * line below, or use the -pc_factor_shift_positive_definite option. 82*c4762a1bSJed Brown * Run the code twice: once as given to see the negative pivot and the 83*c4762a1bSJed Brown * divergence behaviour, then comment in the Shift line, or add the 84*c4762a1bSJed Brown * command line option, and see that the pivots are all positive and 85*c4762a1bSJed Brown * the method converges. 86*c4762a1bSJed Brown */ 87*c4762a1bSJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = PCSetType(pc,PCICC);CHKERRQ(ierr); 89*c4762a1bSJed Brown /* ierr = PCFactorSetShiftType(prec,MAT_SHIFT_POSITIVE_DEFINITE);CHKERRQ(ierr); */ 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = KSPSetUp(ksp);CHKERRQ(ierr); 93*c4762a1bSJed Brown 94*c4762a1bSJed Brown /* 95*c4762a1bSJed Brown * Now that the factorisation is done, show the pivots; 96*c4762a1bSJed Brown * note that the last one is negative. This in itself is not an error, 97*c4762a1bSJed Brown * but it will make the iterative method diverge. 98*c4762a1bSJed Brown */ 99*c4762a1bSJed Brown ierr = PCFactorGetMatrix(pc,&M);CHKERRQ(ierr); 100*c4762a1bSJed Brown ierr = VecDuplicate(B,&D);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = MatGetDiagonal(M,D);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nPivots:\n\n");CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = VecView(D,0);CHKERRQ(ierr); 104*c4762a1bSJed Brown 105*c4762a1bSJed Brown /* 106*c4762a1bSJed Brown * Solve the system; 107*c4762a1bSJed Brown * without the shift this will diverge with 108*c4762a1bSJed Brown * an indefinite preconditioner 109*c4762a1bSJed Brown */ 110*c4762a1bSJed Brown ierr = KSPSolve(ksp,B,X);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = KSPGetConvergedReason(ksp,&reason);CHKERRQ(ierr); 112*c4762a1bSJed Brown if (reason==KSP_DIVERGED_INDEFINITE_PC) { 113*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nDivergence because of indefinite preconditioner;\n");CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Run the executable again but with -pc_factor_shift_positive_definite option.\n");CHKERRQ(ierr); 115*c4762a1bSJed Brown } else if (reason<0) { 116*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nOther kind of divergence: this should not happen.\n");CHKERRQ(ierr); 117*c4762a1bSJed Brown } else { 118*c4762a1bSJed Brown ierr = KSPGetIterationNumber(ksp,&its);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nConvergence in %d iterations.\n",(int)its);CHKERRQ(ierr); 120*c4762a1bSJed Brown } 121*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n");CHKERRQ(ierr); 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 124*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = VecDestroy(&B);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = VecDestroy(&D);CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = PetscFinalize(); 129*c4762a1bSJed Brown return ierr; 130*c4762a1bSJed Brown } 131*c4762a1bSJed Brown 132*c4762a1bSJed Brown 133*c4762a1bSJed Brown /*TEST 134*c4762a1bSJed Brown 135*c4762a1bSJed Brown test: 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown TEST*/ 138