xref: /petsc/src/snes/tutorials/ex59.c (revision 11486bccf1aeea1ca5536228f99d437b39bdaca6) !
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   PetscInt          it,n = 11,i;
33   PetscReal         h,xp = 0.0;
34   PetscScalar       v;
35   const PetscScalar a = X0DOT;
36   const PetscScalar b = X1;
37   const PetscScalar k = KPOW;
38   PetscScalar       v2;
39   PetscScalar       *xx;
40 
41   PetscFunctionBeginUser;
42   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
43   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
44   PetscCall(PetscOptionsGetBool(NULL,NULL,"-second_order",&second_order,NULL));
45   h    = 1.0/(n-1);
46 
47   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48      Create nonlinear solver context
49      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50 
51   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
52 
53   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54      Create vector data structures; set function evaluation routine
55      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56 
57   PetscCall(VecCreate(PETSC_COMM_SELF,&x));
58   PetscCall(VecSetSizes(x,PETSC_DECIDE,n));
59   PetscCall(VecSetFromOptions(x));
60   PetscCall(VecDuplicate(x,&r));
61   PetscCall(VecDuplicate(x,&F));
62 
63   PetscCall(SNESSetFunction(snes,r,FormFunction,(void*)F));
64 
65   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66      Create matrix data structures; set Jacobian evaluation routine
67      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68 
69   PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,n,n,3,NULL,&J));
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   /*  PetscCall(SNESSetJacobian(snes,NULL,JPrec,FormJacobian,0)); */
77   PetscCall(SNESSetJacobian(snes,J,J,FormJacobian,0));
78   /*  PetscCall(SNESSetJacobian(snes,J,JPrec,FormJacobian,0)); */
79 
80   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81      Customize nonlinear solver; set runtime options
82    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83 
84   PetscCall(SNESSetFromOptions(snes));
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     PetscCall(VecSetValues(F,1,&i,&v,INSERT_VALUES));
98     v2   = a*xp + (b-a)*PetscPowScalar(xp,k);
99     PetscCall(VecSetValues(x,1,&i,&v2,INSERT_VALUES));
100     xp  += h;
101   }
102 
103   /* perturb initial guess */
104   PetscCall(VecGetArray(x,&xx));
105   for (i=0; i<n; i++) {
106     v2   = xx[i]*sperturb;
107     PetscCall(VecSetValues(x,1,&i,&v2,INSERT_VALUES));
108   }
109   PetscCall(VecRestoreArray(x,&xx));
110 
111   PetscCall(SNESSolve(snes,NULL,x));
112   PetscCall(SNESGetIterationNumber(snes,&it));
113   PetscCall(PetscPrintf(PETSC_COMM_SELF,"SNES iterations = %" PetscInt_FMT "\n\n",it));
114 
115   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116      Free work space.  All PETSc objects should be destroyed when they
117      are no longer needed.
118    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119 
120   PetscCall(VecDestroy(&x));     PetscCall(VecDestroy(&r));
121   PetscCall(VecDestroy(&F));     PetscCall(MatDestroy(&J));
122   PetscCall(SNESDestroy(&snes));
123   PetscCall(PetscFinalize());
124   return 0;
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   PetscInt          i,n;
132 
133   PetscCall(VecGetArrayRead(x,&xx));
134   PetscCall(VecGetArray(f,&ff));
135   PetscCall(VecGetArray((Vec)dummy,&FF));
136   PetscCall(VecGetSize(x,&n));
137   d    = (PetscReal)(n - 1); d2 = d*d;
138 
139   if (second_order) ff[0] = d*(0.5*d*(-xx[2] + 4.*xx[1] - 3.*xx[0]) - X0DOT);
140   else ff[0] = d*(d*(xx[1] - xx[0]) - X0DOT);
141 
142   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];
143 
144   ff[n-1] = d*d*(xx[n-1] - X1);
145   PetscCall(VecRestoreArrayRead(x,&xx));
146   PetscCall(VecRestoreArray(f,&ff));
147   PetscCall(VecRestoreArray((Vec)dummy,&FF));
148   return 0;
149 }
150 
151 PetscErrorCode FormJacobian(SNES snes,Vec x,Mat jac,Mat prejac,void *dummy)
152 {
153   const PetscScalar *xx;
154   PetscScalar       A[3],d,d2;
155   PetscInt          i,n,j[3];
156 
157   PetscCall(VecGetSize(x,&n));
158   PetscCall(VecGetArrayRead(x,&xx));
159   d    = (PetscReal)(n - 1); d2 = d*d;
160 
161   i = 0;
162   if (second_order) {
163     j[0] = 0; j[1] = 1; j[2] = 2;
164     A[0] = -3.*d*d*0.5; A[1] = 4.*d*d*0.5;  A[2] = -1.*d*d*0.5;
165     PetscCall(MatSetValues(prejac,1,&i,3,j,A,INSERT_VALUES));
166   } else {
167     j[0] = 0; j[1] = 1;
168     A[0] = -d*d; A[1] = d*d;
169     PetscCall(MatSetValues(prejac,1,&i,2,j,A,INSERT_VALUES));
170   }
171   for (i=1; i<n-1; i++) {
172     j[0] = i - 1; j[1] = i;                   j[2] = i + 1;
173     A[0] = d2;    A[1] = -2.*d2 + 2.*xx[i];  A[2] = d2;
174     PetscCall(MatSetValues(prejac,1,&i,3,j,A,INSERT_VALUES));
175   }
176 
177   i    = n-1;
178   A[0] = d*d;
179   PetscCall(MatSetValues(prejac,1,&i,1,&i,&A[0],INSERT_VALUES));
180 
181   PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
182   PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
183   PetscCall(MatAssemblyBegin(prejac,MAT_FINAL_ASSEMBLY));
184   PetscCall(MatAssemblyEnd(prejac,MAT_FINAL_ASSEMBLY));
185 
186   PetscCall(VecRestoreArrayRead(x,&xx));
187   return 0;
188 }
189 
190 /*TEST
191 
192    test:
193       args: -n 14 -snes_monitor_short -snes_converged_reason
194       requires: !single
195 
196    test:
197       suffix: 2
198       args: -n 15 -snes_monitor_short -snes_converged_reason
199       requires: !single
200 
201    test:
202       suffix: 3
203       args: -n 14 -second_order -snes_monitor_short -snes_converged_reason
204       requires: !single
205 
206 TEST*/
207