1 2 static const char help[] = "Tries to solve u`` + u^{2} = f for an easy case and an impossible case.\n\n"; 3 4 /* 5 This example was contributed by Peter Graf to show how SNES fails when given a nonlinear problem with no solution. 6 7 Run with -n 14 to see fail to converge and -n 15 to see convergence 8 9 The option -second_order uses a different discretization of the Neumann boundary condition and always converges 10 11 */ 12 13 #include <petscsnes.h> 14 15 PetscBool second_order = PETSC_FALSE; 16 #define X0DOT -2.0 17 #define X1 5.0 18 #define KPOW 2.0 19 const PetscScalar sperturb = 1.1; 20 21 /* 22 User-defined routines 23 */ 24 PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*); 25 PetscErrorCode FormFunction(SNES,Vec,Vec,void*); 26 27 int main(int argc,char **argv) 28 { 29 SNES snes; /* SNES context */ 30 Vec x,r,F; /* vectors */ 31 Mat J; /* Jacobian */ 32 PetscErrorCode ierr; 33 PetscInt it,n = 11,i; 34 PetscReal h,xp = 0.0; 35 PetscScalar v; 36 const PetscScalar a = X0DOT; 37 const PetscScalar b = X1; 38 const PetscScalar k = KPOW; 39 PetscScalar v2; 40 PetscScalar *xx; 41 42 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 43 ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 44 ierr = PetscOptionsGetBool(NULL,NULL,"-second_order",&second_order,NULL);CHKERRQ(ierr); 45 h = 1.0/(n-1); 46 47 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 48 Create nonlinear solver context 49 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 50 51 ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 52 53 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 54 Create vector data structures; set function evaluation routine 55 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 56 57 ierr = VecCreate(PETSC_COMM_SELF,&x);CHKERRQ(ierr); 58 ierr = VecSetSizes(x,PETSC_DECIDE,n);CHKERRQ(ierr); 59 ierr = VecSetFromOptions(x);CHKERRQ(ierr); 60 ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 61 ierr = VecDuplicate(x,&F);CHKERRQ(ierr); 62 63 ierr = SNESSetFunction(snes,r,FormFunction,(void*)F);CHKERRQ(ierr); 64 65 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 66 Create matrix data structures; set Jacobian evaluation routine 67 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 68 69 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&J);CHKERRQ(ierr); 70 71 /* 72 Note that in this case we create separate matrices for the Jacobian 73 and preconditioner matrix. Both of these are computed in the 74 routine FormJacobian() 75 */ 76 /* ierr = SNESSetJacobian(snes,NULL,JPrec,FormJacobian,0);CHKERRQ(ierr); */ 77 ierr = SNESSetJacobian(snes,J,J,FormJacobian,0);CHKERRQ(ierr); 78 /* ierr = SNESSetJacobian(snes,J,JPrec,FormJacobian,0);CHKERRQ(ierr); */ 79 80 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 81 Customize nonlinear solver; set runtime options 82 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 83 84 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 85 86 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 87 Initialize application: 88 Store right-hand-side of PDE and exact solution 89 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 90 91 /* set right hand side and initial guess to be exact solution of continuim problem */ 92 #define SQR(x) ((x)*(x)) 93 xp = 0.0; 94 for (i=0; i<n; i++) 95 { 96 v = k*(k-1.)*(b-a)*PetscPowScalar(xp,k-2.) + SQR(a*xp) + SQR(b-a)*PetscPowScalar(xp,2.*k) + 2.*a*(b-a)*PetscPowScalar(xp,k+1.); 97 ierr = VecSetValues(F,1,&i,&v,INSERT_VALUES);CHKERRQ(ierr); 98 v2 = a*xp + (b-a)*PetscPowScalar(xp,k); 99 ierr = VecSetValues(x,1,&i,&v2,INSERT_VALUES);CHKERRQ(ierr); 100 xp += h; 101 } 102 103 /* perturb initial guess */ 104 ierr = VecGetArray(x,&xx);CHKERRQ(ierr); 105 for (i=0; i<n; i++) { 106 v2 = xx[i]*sperturb; 107 ierr = VecSetValues(x,1,&i,&v2,INSERT_VALUES);CHKERRQ(ierr); 108 } 109 ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr); 110 111 ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr); 112 ierr = SNESGetIterationNumber(snes,&it);CHKERRQ(ierr); 113 ierr = PetscPrintf(PETSC_COMM_SELF,"SNES iterations = %D\n\n",it);CHKERRQ(ierr); 114 115 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116 Free work space. All PETSc objects should be destroyed when they 117 are no longer needed. 118 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 119 120 ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr); 121 ierr = VecDestroy(&F);CHKERRQ(ierr); ierr = MatDestroy(&J);CHKERRQ(ierr); 122 ierr = SNESDestroy(&snes);CHKERRQ(ierr); 123 ierr = PetscFinalize(); 124 return ierr; 125 } 126 127 PetscErrorCode FormFunction(SNES snes,Vec x,Vec f,void *dummy) 128 { 129 const PetscScalar *xx; 130 PetscScalar *ff,*FF,d,d2; 131 PetscErrorCode ierr; 132 PetscInt i,n; 133 134 ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 135 ierr = VecGetArray(f,&ff);CHKERRQ(ierr); 136 ierr = VecGetArray((Vec)dummy,&FF);CHKERRQ(ierr); 137 ierr = VecGetSize(x,&n);CHKERRQ(ierr); 138 d = (PetscReal)(n - 1); d2 = d*d; 139 140 if (second_order) ff[0] = d*(0.5*d*(-xx[2] + 4.*xx[1] - 3.*xx[0]) - X0DOT); 141 else ff[0] = d*(d*(xx[1] - xx[0]) - X0DOT); 142 143 for (i=1; i<n-1; i++) ff[i] = d2*(xx[i-1] - 2.*xx[i] + xx[i+1]) + xx[i]*xx[i] - FF[i]; 144 145 ff[n-1] = d*d*(xx[n-1] - X1); 146 ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 147 ierr = VecRestoreArray(f,&ff);CHKERRQ(ierr); 148 ierr = VecRestoreArray((Vec)dummy,&FF);CHKERRQ(ierr); 149 return 0; 150 } 151 152 PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat prejac,void *dummy) 153 { 154 const PetscScalar *xx; 155 PetscScalar A[3],d,d2; 156 PetscInt i,n,j[3]; 157 PetscErrorCode ierr; 158 159 ierr = VecGetSize(x,&n);CHKERRQ(ierr); 160 ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 161 d = (PetscReal)(n - 1); d2 = d*d; 162 163 i = 0; 164 if (second_order) { 165 j[0] = 0; j[1] = 1; j[2] = 2; 166 A[0] = -3.*d*d*0.5; A[1] = 4.*d*d*0.5; A[2] = -1.*d*d*0.5; 167 ierr = MatSetValues(prejac,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr); 168 } else { 169 j[0] = 0; j[1] = 1; 170 A[0] = -d*d; A[1] = d*d; 171 ierr = MatSetValues(prejac,1,&i,2,j,A,INSERT_VALUES);CHKERRQ(ierr); 172 } 173 for (i=1; i<n-1; i++) { 174 j[0] = i - 1; j[1] = i; j[2] = i + 1; 175 A[0] = d2; A[1] = -2.*d2 + 2.*xx[i]; A[2] = d2; 176 ierr = MatSetValues(prejac,1,&i,3,j,A,INSERT_VALUES);CHKERRQ(ierr); 177 } 178 179 i = n-1; 180 A[0] = d*d; 181 ierr = MatSetValues(prejac,1,&i,1,&i,&A[0],INSERT_VALUES);CHKERRQ(ierr); 182 183 ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 184 ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 185 ierr = MatAssemblyBegin(prejac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 186 ierr = MatAssemblyEnd(prejac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 187 188 ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 189 return 0; 190 } 191 192 /*TEST 193 194 test: 195 args: -n 14 -snes_monitor_short -snes_converged_reason 196 requires: !single 197 198 test: 199 suffix: 2 200 args: -n 15 -snes_monitor_short -snes_converged_reason 201 requires: !single 202 203 test: 204 suffix: 3 205 args: -n 14 -second_order -snes_monitor_short -snes_converged_reason 206 requires: !single 207 208 TEST*/ 209