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_type 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 * Contributed by Victor Eijkhout 2003. 16*c4762a1bSJed Brown */ 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown #include <petscksp.h> 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown int main(int argc,char **argv) 21*c4762a1bSJed Brown { 22*c4762a1bSJed Brown KSP solver; 23*c4762a1bSJed Brown PC prec; 24*c4762a1bSJed Brown Mat A,M; 25*c4762a1bSJed Brown Vec X,B,D; 26*c4762a1bSJed Brown MPI_Comm comm; 27*c4762a1bSJed Brown PetscScalar v; 28*c4762a1bSJed Brown KSPConvergedReason reason; 29*c4762a1bSJed Brown PetscInt i,j,its; 30*c4762a1bSJed Brown PetscErrorCode ierr; 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 33*c4762a1bSJed Brown comm = MPI_COMM_SELF; 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown /* 36*c4762a1bSJed Brown * Construct the Kershaw matrix 37*c4762a1bSJed Brown * and a suitable rhs / initial guess 38*c4762a1bSJed Brown */ 39*c4762a1bSJed Brown ierr = MatCreateSeqAIJ(comm,4,4,4,0,&A);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = VecCreateSeq(comm,4,&B);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = VecDuplicate(B,&X);CHKERRQ(ierr); 42*c4762a1bSJed Brown for (i=0; i<4; i++) { 43*c4762a1bSJed Brown v = 3; 44*c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 45*c4762a1bSJed Brown v = 1; 46*c4762a1bSJed Brown ierr = VecSetValues(B,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = VecSetValues(X,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 48*c4762a1bSJed Brown } 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown i=0; v=0; 51*c4762a1bSJed Brown ierr = VecSetValues(X,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown for (i=0; i<3; i++) { 54*c4762a1bSJed Brown v = -2; j=i+1; 55*c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 57*c4762a1bSJed Brown } 58*c4762a1bSJed Brown i=0; j=3; v=2; 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = MatSetValues(A,1,&j,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = VecAssemblyBegin(B);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = VecAssemblyEnd(B);CHKERRQ(ierr); 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown /* 68*c4762a1bSJed Brown * A Conjugate Gradient method 69*c4762a1bSJed Brown * with ILU(0) preconditioning 70*c4762a1bSJed Brown */ 71*c4762a1bSJed Brown ierr = KSPCreate(comm,&solver);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = KSPSetOperators(solver,A,A);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown ierr = KSPSetType(solver,KSPCG);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = KSPSetInitialGuessNonzero(solver,PETSC_TRUE);CHKERRQ(ierr); 76*c4762a1bSJed Brown 77*c4762a1bSJed Brown /* 78*c4762a1bSJed Brown * ILU preconditioner; 79*c4762a1bSJed Brown * this will break down unless you add the Shift line, 80*c4762a1bSJed Brown * or use the -pc_factor_shift_positive_definite option */ 81*c4762a1bSJed Brown ierr = KSPGetPC(solver,&prec);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = PCSetType(prec,PCILU);CHKERRQ(ierr); 83*c4762a1bSJed Brown /* ierr = PCFactorSetShiftType(prec,MAT_SHIFT_POSITIVE_DEFINITE);CHKERRQ(ierr); */ 84*c4762a1bSJed Brown 85*c4762a1bSJed Brown ierr = KSPSetFromOptions(solver);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = KSPSetUp(solver);CHKERRQ(ierr); 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown /* 89*c4762a1bSJed Brown * Now that the factorisation is done, show the pivots; 90*c4762a1bSJed Brown * note that the last one is negative. This in itself is not an error, 91*c4762a1bSJed Brown * but it will make the iterative method diverge. 92*c4762a1bSJed Brown */ 93*c4762a1bSJed Brown ierr = PCFactorGetMatrix(prec,&M);CHKERRQ(ierr); 94*c4762a1bSJed Brown ierr = VecDuplicate(B,&D);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = MatGetDiagonal(M,D);CHKERRQ(ierr); 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown /* 98*c4762a1bSJed Brown * Solve the system; 99*c4762a1bSJed Brown * without the shift this will diverge with 100*c4762a1bSJed Brown * an indefinite preconditioner 101*c4762a1bSJed Brown */ 102*c4762a1bSJed Brown ierr = KSPSolve(solver,B,X);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = KSPGetConvergedReason(solver,&reason);CHKERRQ(ierr); 104*c4762a1bSJed Brown if (reason==KSP_DIVERGED_INDEFINITE_PC) { 105*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nDivergence because of indefinite preconditioner;\n");CHKERRQ(ierr); 106*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Run the executable again but with '-pc_factor_shift_type POSITIVE_DEFINITE' option.\n");CHKERRQ(ierr); 107*c4762a1bSJed Brown } else if (reason<0) { 108*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nOther kind of divergence: this should not happen.\n");CHKERRQ(ierr); 109*c4762a1bSJed Brown } else { 110*c4762a1bSJed Brown ierr = KSPGetIterationNumber(solver,&its);CHKERRQ(ierr); 111*c4762a1bSJed Brown } 112*c4762a1bSJed Brown 113*c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = VecDestroy(&B);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = VecDestroy(&D);CHKERRQ(ierr); 116*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 117*c4762a1bSJed Brown ierr = KSPDestroy(&solver);CHKERRQ(ierr); 118*c4762a1bSJed Brown ierr = PetscFinalize(); 119*c4762a1bSJed Brown return ierr; 120*c4762a1bSJed Brown } 121*c4762a1bSJed Brown 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown /*TEST 124*c4762a1bSJed Brown 125*c4762a1bSJed Brown test: 126*c4762a1bSJed Brown args: -pc_factor_shift_type positive_definite 127*c4762a1bSJed Brown 128*c4762a1bSJed Brown TEST*/ 129