xref: /petsc/src/snes/tutorials/ex55.c (revision f5a72eb34e1dc21bf8593cf383cc1fa498006d28)
1 static char help[] = "Copy of ex5.c\n";
2 
3 /* ------------------------------------------------------------------------
4 
5   Copy of ex5.c.
6   Once petsc test harness supports conditional linking, we can remove this duplicate.
7   See https://gitlab.com/petsc/petsc/-/issues/1173
8   ------------------------------------------------------------------------- */
9 
10 /*
11    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
12    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
13 */
14 #include <petscdm.h>
15 #include <petscdmda.h>
16 #include <petscsnes.h>
17 #include <petscmatlab.h>
18 #include <petsc/private/snesimpl.h> /* For SNES_Solve event */
19 #include "ex55.h"
20 
21 /* ------------------------------------------------------------------- */
22 /*
23    FormInitialGuess - Forms initial approximation.
24 
25    Input Parameters:
26    da - The DM
27    user - user-defined application context
28 
29    Output Parameter:
30    X - vector
31  */
32 static PetscErrorCode FormInitialGuess(DM da,AppCtx *user,Vec X)
33 {
34   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
35   PetscReal      lambda,temp1,temp,hx,hy;
36   PetscScalar    **x;
37 
38   PetscFunctionBeginUser;
39   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
40 
41   lambda = user->param;
42   hx     = 1.0/(PetscReal)(Mx-1);
43   hy     = 1.0/(PetscReal)(My-1);
44   temp1  = lambda/(lambda + 1.0);
45 
46   /*
47      Get a pointer to vector data.
48        - For default PETSc vectors, VecGetArray() returns a pointer to
49          the data array.  Otherwise, the routine is implementation dependent.
50        - You MUST call VecRestoreArray() when you no longer need access to
51          the array.
52   */
53   PetscCall(DMDAVecGetArray(da,X,&x));
54 
55   /*
56      Get local grid boundaries (for 2-dimensional DMDA):
57        xs, ys   - starting grid indices (no ghost points)
58        xm, ym   - widths of local grid (no ghost points)
59 
60   */
61   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
62 
63   /*
64      Compute initial guess over the locally owned part of the grid
65   */
66   for (j=ys; j<ys+ym; j++) {
67     temp = (PetscReal)(PetscMin(j,My-j-1))*hy;
68     for (i=xs; i<xs+xm; i++) {
69       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
70         /* boundary conditions are all zero Dirichlet */
71         x[j][i] = 0.0;
72       } else {
73         x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
74       }
75     }
76   }
77 
78   /*
79      Restore vector
80   */
81   PetscCall(DMDAVecRestoreArray(da,X,&x));
82   PetscFunctionReturn(0);
83 }
84 
85 /*
86   FormExactSolution - Forms MMS solution
87 
88   Input Parameters:
89   da - The DM
90   user - user-defined application context
91 
92   Output Parameter:
93   X - vector
94  */
95 static PetscErrorCode FormExactSolution(DM da, AppCtx *user, Vec U)
96 {
97   DM             coordDA;
98   Vec            coordinates;
99   DMDACoor2d   **coords;
100   PetscScalar  **u;
101   PetscInt       xs, ys, xm, ym, i, j;
102 
103   PetscFunctionBeginUser;
104   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
105   PetscCall(DMGetCoordinateDM(da, &coordDA));
106   PetscCall(DMGetCoordinates(da, &coordinates));
107   PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords));
108   PetscCall(DMDAVecGetArray(da, U, &u));
109   for (j = ys; j < ys+ym; ++j) {
110     for (i = xs; i < xs+xm; ++i) {
111       user->mms_solution(user,&coords[j][i],&u[j][i]);
112     }
113   }
114   PetscCall(DMDAVecRestoreArray(da, U, &u));
115   PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords));
116   PetscFunctionReturn(0);
117 }
118 
119 static PetscErrorCode ZeroBCSolution(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
120 {
121   u[0] = 0.;
122   return 0;
123 }
124 
125 /* The functions below evaluate the MMS solution u(x,y) and associated forcing
126 
127      f(x,y) = -u_xx - u_yy - lambda exp(u)
128 
129   such that u(x,y) is an exact solution with f(x,y) as the right hand side forcing term.
130  */
131 static PetscErrorCode MMSSolution1(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
132 {
133   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
134   u[0] = x*(1 - x)*y*(1 - y);
135   PetscLogFlops(5);
136   return 0;
137 }
138 static PetscErrorCode MMSForcing1(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
139 {
140   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
141   f[0] = 2*x*(1 - x) + 2*y*(1 - y) - user->param*PetscExpReal(x*(1 - x)*y*(1 - y));
142   return 0;
143 }
144 
145 static PetscErrorCode MMSSolution2(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
146 {
147   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
148   u[0] = PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y);
149   PetscLogFlops(5);
150   return 0;
151 }
152 static PetscErrorCode MMSForcing2(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
153 {
154   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
155   f[0] = 2*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y) - user->param*PetscExpReal(PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y));
156   return 0;
157 }
158 
159 static PetscErrorCode MMSSolution3(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
160 {
161   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
162   u[0] = PetscSinReal(user->m*PETSC_PI*x*(1-y))*PetscSinReal(user->n*PETSC_PI*y*(1-x));
163   PetscLogFlops(5);
164   return 0;
165 }
166 static PetscErrorCode MMSForcing3(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
167 {
168   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
169   PetscReal m = user->m, n = user->n, lambda = user->param;
170   f[0] = (-(PetscExpReal(PetscSinReal(m*PETSC_PI*x*(1 - y))*PetscSinReal(n*PETSC_PI*(1 - x)*y))*lambda)
171           + PetscSqr(PETSC_PI)*(-2*m*n*((-1 + x)*x + (-1 + y)*y)*PetscCosReal(m*PETSC_PI*x*(-1 + y))*PetscCosReal(n*PETSC_PI*(-1 + x)*y)
172                                 + (PetscSqr(m)*(PetscSqr(x) + PetscSqr(-1 + y)) + PetscSqr(n)*(PetscSqr(-1 + x) + PetscSqr(y)))
173                                 *PetscSinReal(m*PETSC_PI*x*(-1 + y))*PetscSinReal(n*PETSC_PI*(-1 + x)*y)));
174   return 0;
175 }
176 
177 static PetscErrorCode MMSSolution4(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
178 {
179   const PetscReal Lx = 1.,Ly = 1.;
180   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
181   u[0] = (PetscPowReal(x,4)-PetscSqr(Lx)*PetscSqr(x))*(PetscPowReal(y,4)-PetscSqr(Ly)*PetscSqr(y));
182   PetscLogFlops(9);
183   return 0;
184 }
185 static PetscErrorCode MMSForcing4(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
186 {
187   const PetscReal Lx = 1.,Ly = 1.;
188   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
189   f[0] = (2*PetscSqr(x)*(PetscSqr(x)-PetscSqr(Lx))*(PetscSqr(Ly)-6*PetscSqr(y))
190           + 2*PetscSqr(y)*(PetscSqr(Lx)-6*PetscSqr(x))*(PetscSqr(y)-PetscSqr(Ly))
191           - user->param*PetscExpReal((PetscPowReal(x,4)-PetscSqr(Lx)*PetscSqr(x))*(PetscPowReal(y,4)-PetscSqr(Ly)*PetscSqr(y))));
192   return 0;
193 }
194 
195 /* ------------------------------------------------------------------- */
196 /*
197    FormFunctionLocal - Evaluates nonlinear function, F(x) on local process patch
198 
199  */
200 static PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
201 {
202   PetscInt       i,j;
203   PetscReal      lambda,hx,hy,hxdhy,hydhx;
204   PetscScalar    u,ue,uw,un,us,uxx,uyy,mms_solution,mms_forcing;
205   DMDACoor2d     c;
206 
207   PetscFunctionBeginUser;
208   lambda = user->param;
209   hx     = 1.0/(PetscReal)(info->mx-1);
210   hy     = 1.0/(PetscReal)(info->my-1);
211   hxdhy  = hx/hy;
212   hydhx  = hy/hx;
213   /*
214      Compute function over the locally owned part of the grid
215   */
216   for (j=info->ys; j<info->ys+info->ym; j++) {
217     for (i=info->xs; i<info->xs+info->xm; i++) {
218       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
219         c.x = i*hx; c.y = j*hy;
220         PetscCall(user->mms_solution(user,&c,&mms_solution));
221         f[j][i] = 2.0*(hydhx+hxdhy)*(x[j][i] - mms_solution);
222       } else {
223         u  = x[j][i];
224         uw = x[j][i-1];
225         ue = x[j][i+1];
226         un = x[j-1][i];
227         us = x[j+1][i];
228 
229         /* Enforce boundary conditions at neighboring points -- setting these values causes the Jacobian to be symmetric. */
230         if (i-1 == 0) {c.x = (i-1)*hx; c.y = j*hy; PetscCall(user->mms_solution(user,&c,&uw));}
231         if (i+1 == info->mx-1) {c.x = (i+1)*hx; c.y = j*hy; PetscCall(user->mms_solution(user,&c,&ue));}
232         if (j-1 == 0) {c.x = i*hx; c.y = (j-1)*hy; PetscCall(user->mms_solution(user,&c,&un));}
233         if (j+1 == info->my-1) {c.x = i*hx; c.y = (j+1)*hy; PetscCall(user->mms_solution(user,&c,&us));}
234 
235         uxx     = (2.0*u - uw - ue)*hydhx;
236         uyy     = (2.0*u - un - us)*hxdhy;
237         mms_forcing = 0;
238         c.x = i*hx; c.y = j*hy;
239         if (user->mms_forcing) PetscCall(user->mms_forcing(user,&c,&mms_forcing));
240         f[j][i] = uxx + uyy - hx*hy*(lambda*PetscExpScalar(u) + mms_forcing);
241       }
242     }
243   }
244   PetscCall(PetscLogFlops(11.0*info->ym*info->xm));
245   PetscFunctionReturn(0);
246 }
247 
248 /* FormObjectiveLocal - Evaluates nonlinear function, F(x) on local process patch */
249 static PetscErrorCode FormObjectiveLocal(DMDALocalInfo *info,PetscScalar **x,PetscReal *obj,AppCtx *user)
250 {
251   PetscInt       i,j;
252   PetscReal      lambda,hx,hy,hxdhy,hydhx,sc,lobj=0;
253   PetscScalar    u,ue,uw,un,us,uxux,uyuy;
254   MPI_Comm       comm;
255 
256   PetscFunctionBeginUser;
257   *obj   = 0;
258   PetscCall(PetscObjectGetComm((PetscObject)info->da,&comm));
259   lambda = user->param;
260   hx     = 1.0/(PetscReal)(info->mx-1);
261   hy     = 1.0/(PetscReal)(info->my-1);
262   sc     = hx*hy*lambda;
263   hxdhy  = hx/hy;
264   hydhx  = hy/hx;
265   /*
266      Compute function over the locally owned part of the grid
267   */
268   for (j=info->ys; j<info->ys+info->ym; j++) {
269     for (i=info->xs; i<info->xs+info->xm; i++) {
270       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
271         lobj += PetscRealPart((hydhx + hxdhy)*x[j][i]*x[j][i]);
272       } else {
273         u  = x[j][i];
274         uw = x[j][i-1];
275         ue = x[j][i+1];
276         un = x[j-1][i];
277         us = x[j+1][i];
278 
279         if (i-1 == 0) uw = 0.;
280         if (i+1 == info->mx-1) ue = 0.;
281         if (j-1 == 0) un = 0.;
282         if (j+1 == info->my-1) us = 0.;
283 
284         /* F[u] = 1/2\int_{\omega}\nabla^2u(x)*u(x)*dx */
285 
286         uxux = u*(2.*u - ue - uw)*hydhx;
287         uyuy = u*(2.*u - un - us)*hxdhy;
288 
289         lobj += PetscRealPart(0.5*(uxux + uyuy) - sc*PetscExpScalar(u));
290       }
291     }
292   }
293   PetscCall(PetscLogFlops(12.0*info->ym*info->xm));
294   PetscCallMPI(MPI_Allreduce(&lobj,obj,1,MPIU_REAL,MPIU_SUM,comm));
295   PetscFunctionReturn(0);
296 }
297 
298 /*
299    FormJacobianLocal - Evaluates Jacobian matrix on local process patch
300 */
301 static PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat jac,Mat jacpre,AppCtx *user)
302 {
303   PetscInt       i,j,k;
304   MatStencil     col[5],row;
305   PetscScalar    lambda,v[5],hx,hy,hxdhy,hydhx,sc;
306   DM             coordDA;
307   Vec            coordinates;
308   DMDACoor2d   **coords;
309 
310   PetscFunctionBeginUser;
311   lambda = user->param;
312   /* Extract coordinates */
313   PetscCall(DMGetCoordinateDM(info->da, &coordDA));
314   PetscCall(DMGetCoordinates(info->da, &coordinates));
315   PetscCall(DMDAVecGetArray(coordDA, coordinates, &coords));
316   hx     = info->xm > 1 ? PetscRealPart(coords[info->ys][info->xs+1].x) - PetscRealPart(coords[info->ys][info->xs].x) : 1.0;
317   hy     = info->ym > 1 ? PetscRealPart(coords[info->ys+1][info->xs].y) - PetscRealPart(coords[info->ys][info->xs].y) : 1.0;
318   PetscCall(DMDAVecRestoreArray(coordDA, coordinates, &coords));
319   hxdhy  = hx/hy;
320   hydhx  = hy/hx;
321   sc     = hx*hy*lambda;
322 
323   /*
324      Compute entries for the locally owned part of the Jacobian.
325       - Currently, all PETSc parallel matrix formats are partitioned by
326         contiguous chunks of rows across the processors.
327       - Each processor needs to insert only elements that it owns
328         locally (but any non-local elements will be sent to the
329         appropriate processor during matrix assembly).
330       - Here, we set all entries for a particular row at once.
331       - We can set matrix entries either using either
332         MatSetValuesLocal() or MatSetValues(), as discussed above.
333   */
334   for (j=info->ys; j<info->ys+info->ym; j++) {
335     for (i=info->xs; i<info->xs+info->xm; i++) {
336       row.j = j; row.i = i;
337       /* boundary points */
338       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
339         v[0] =  2.0*(hydhx + hxdhy);
340         PetscCall(MatSetValuesStencil(jacpre,1,&row,1,&row,v,INSERT_VALUES));
341       } else {
342         k = 0;
343         /* interior grid points */
344         if (j-1 != 0) {
345           v[k]     = -hxdhy;
346           col[k].j = j - 1; col[k].i = i;
347           k++;
348         }
349         if (i-1 != 0) {
350           v[k]     = -hydhx;
351           col[k].j = j;     col[k].i = i-1;
352           k++;
353         }
354 
355         v[k] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(x[j][i]); col[k].j = row.j; col[k].i = row.i; k++;
356 
357         if (i+1 != info->mx-1) {
358           v[k]     = -hydhx;
359           col[k].j = j;     col[k].i = i+1;
360           k++;
361         }
362         if (j+1 != info->mx-1) {
363           v[k]     = -hxdhy;
364           col[k].j = j + 1; col[k].i = i;
365           k++;
366         }
367         PetscCall(MatSetValuesStencil(jacpre,1,&row,k,col,v,INSERT_VALUES));
368       }
369     }
370   }
371 
372   /*
373      Assemble matrix, using the 2-step process:
374        MatAssemblyBegin(), MatAssemblyEnd().
375   */
376   PetscCall(MatAssemblyBegin(jacpre,MAT_FINAL_ASSEMBLY));
377   PetscCall(MatAssemblyEnd(jacpre,MAT_FINAL_ASSEMBLY));
378 
379   /*
380      Tell the matrix we will never add a new nonzero location to the
381      matrix. If we do, it will generate an error.
382   */
383   PetscCall(MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
384   PetscFunctionReturn(0);
385 }
386 
387 static PetscErrorCode FormFunctionMatlab(SNES snes,Vec X,Vec F,void *ptr)
388 {
389 #if PetscDefined(HAVE_MATLAB_ENGINE)
390   AppCtx         *user = (AppCtx*)ptr;
391   PetscInt       Mx,My;
392   PetscReal      lambda,hx,hy;
393   Vec            localX,localF;
394   MPI_Comm       comm;
395   DM             da;
396 
397   PetscFunctionBeginUser;
398   PetscCall(SNESGetDM(snes,&da));
399   PetscCall(DMGetLocalVector(da,&localX));
400   PetscCall(DMGetLocalVector(da,&localF));
401   PetscCall(PetscObjectSetName((PetscObject)localX,"localX"));
402   PetscCall(PetscObjectSetName((PetscObject)localF,"localF"));
403   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
404 
405   lambda = user->param;
406   hx     = 1.0/(PetscReal)(Mx-1);
407   hy     = 1.0/(PetscReal)(My-1);
408 
409   PetscCall(PetscObjectGetComm((PetscObject)snes,&comm));
410   /*
411      Scatter ghost points to local vector,using the 2-step process
412         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
413      By placing code between these two statements, computations can be
414      done while messages are in transition.
415   */
416   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
417   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
418   PetscCall(PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localX));
419   PetscCall(PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm),"localF=ex5m(localX,%18.16e,%18.16e,%18.16e)",(double)hx,(double)hy,(double)lambda));
420   PetscCall(PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localF));
421 
422   /*
423      Insert values into global vector
424   */
425   PetscCall(DMLocalToGlobalBegin(da,localF,INSERT_VALUES,F));
426   PetscCall(DMLocalToGlobalEnd(da,localF,INSERT_VALUES,F));
427   PetscCall(DMRestoreLocalVector(da,&localX));
428   PetscCall(DMRestoreLocalVector(da,&localF));
429   PetscFunctionReturn(0);
430 #else
431     return 0;                     /* Never called */
432 #endif
433 }
434 
435 /* ------------------------------------------------------------------- */
436 /*
437       Applies some sweeps on nonlinear Gauss-Seidel on each process
438 
439  */
440 static PetscErrorCode NonlinearGS(SNES snes,Vec X, Vec B, void *ctx)
441 {
442   PetscInt       i,j,k,Mx,My,xs,ys,xm,ym,its,tot_its,sweeps,l;
443   PetscReal      lambda,hx,hy,hxdhy,hydhx,sc;
444   PetscScalar    **x,**b,bij,F,F0=0,J,u,un,us,ue,eu,uw,uxx,uyy,y;
445   PetscReal      atol,rtol,stol;
446   DM             da;
447   AppCtx         *user;
448   Vec            localX,localB;
449 
450   PetscFunctionBeginUser;
451   tot_its = 0;
452   PetscCall(SNESNGSGetSweeps(snes,&sweeps));
453   PetscCall(SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its));
454   PetscCall(SNESGetDM(snes,&da));
455   PetscCall(DMGetApplicationContext(da,&user));
456 
457   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
458 
459   lambda = user->param;
460   hx     = 1.0/(PetscReal)(Mx-1);
461   hy     = 1.0/(PetscReal)(My-1);
462   sc     = hx*hy*lambda;
463   hxdhy  = hx/hy;
464   hydhx  = hy/hx;
465 
466   PetscCall(DMGetLocalVector(da,&localX));
467   if (B) {
468     PetscCall(DMGetLocalVector(da,&localB));
469   }
470   for (l=0; l<sweeps; l++) {
471 
472     PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
473     PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
474     if (B) {
475       PetscCall(DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB));
476       PetscCall(DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB));
477     }
478     /*
479      Get a pointer to vector data.
480      - For default PETSc vectors, VecGetArray() returns a pointer to
481      the data array.  Otherwise, the routine is implementation dependent.
482      - You MUST call VecRestoreArray() when you no longer need access to
483      the array.
484      */
485     PetscCall(DMDAVecGetArray(da,localX,&x));
486     if (B) PetscCall(DMDAVecGetArray(da,localB,&b));
487     /*
488      Get local grid boundaries (for 2-dimensional DMDA):
489      xs, ys   - starting grid indices (no ghost points)
490      xm, ym   - widths of local grid (no ghost points)
491      */
492     PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
493 
494     for (j=ys; j<ys+ym; j++) {
495       for (i=xs; i<xs+xm; i++) {
496         if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
497           /* boundary conditions are all zero Dirichlet */
498           x[j][i] = 0.0;
499         } else {
500           if (B) bij = b[j][i];
501           else   bij = 0.;
502 
503           u  = x[j][i];
504           un = x[j-1][i];
505           us = x[j+1][i];
506           ue = x[j][i-1];
507           uw = x[j][i+1];
508 
509           for (k=0; k<its; k++) {
510             eu  = PetscExpScalar(u);
511             uxx = (2.0*u - ue - uw)*hydhx;
512             uyy = (2.0*u - un - us)*hxdhy;
513             F   = uxx + uyy - sc*eu - bij;
514             if (k == 0) F0 = F;
515             J  = 2.0*(hydhx + hxdhy) - sc*eu;
516             y  = F/J;
517             u -= y;
518             tot_its++;
519 
520             if (atol > PetscAbsReal(PetscRealPart(F)) ||
521                 rtol*PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) ||
522                 stol*PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) {
523               break;
524             }
525           }
526           x[j][i] = u;
527         }
528       }
529     }
530     /*
531      Restore vector
532      */
533     PetscCall(DMDAVecRestoreArray(da,localX,&x));
534     PetscCall(DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X));
535     PetscCall(DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X));
536   }
537   PetscCall(PetscLogFlops(tot_its*(21.0)));
538   PetscCall(DMRestoreLocalVector(da,&localX));
539   if (B) {
540     PetscCall(DMDAVecRestoreArray(da,localB,&b));
541     PetscCall(DMRestoreLocalVector(da,&localB));
542   }
543   PetscFunctionReturn(0);
544 }
545 
546 int main(int argc,char **argv)
547 {
548   SNES           snes;                         /* nonlinear solver */
549   Vec            x;                            /* solution vector */
550   AppCtx         user;                         /* user-defined work context */
551   PetscInt       its;                          /* iterations for convergence */
552   PetscReal      bratu_lambda_max = 6.81;
553   PetscReal      bratu_lambda_min = 0.;
554   PetscInt       MMS              = 1;
555   PetscBool      flg              = PETSC_FALSE,setMMS;
556   DM             da;
557   Vec            r                = NULL;
558   KSP            ksp;
559   PetscInt       lits,slits;
560   PetscBool      useKokkos        = PETSC_FALSE;
561 
562   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
563      Initialize program
564      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
565 
566   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
567 
568   PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_kokkos",&useKokkos,NULL));
569 
570   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
571      Initialize problem parameters
572   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
573   user.param = 6.0;
574   PetscCall(PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL));
575   PetscCheck(user.param <= bratu_lambda_max && user.param >= bratu_lambda_min,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lambda, %g, is out of range, [%g, %g]", (double)user.param, (double)bratu_lambda_min, (double)bratu_lambda_max);
576   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mms",&MMS,&setMMS));
577   if (MMS == 3) {
578     PetscInt mPar = 2, nPar = 1;
579     PetscCall(PetscOptionsGetInt(NULL,NULL,"-m_par",&mPar,NULL));
580     PetscCall(PetscOptionsGetInt(NULL,NULL,"-n_par",&nPar,NULL));
581     user.m = PetscPowInt(2,mPar);
582     user.n = PetscPowInt(2,nPar);
583   }
584 
585   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
586      Create nonlinear solver context
587      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
588   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
589   PetscCall(SNESSetCountersReset(snes,PETSC_FALSE));
590   PetscCall(SNESSetNGS(snes, NonlinearGS, NULL));
591 
592   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
593      Create distributed array (DMDA) to manage parallel grid and vectors
594   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
595   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
596   PetscCall(DMSetFromOptions(da));
597   PetscCall(DMSetUp(da));
598   PetscCall(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
599   PetscCall(DMSetApplicationContext(da,&user));
600   PetscCall(SNESSetDM(snes,da));
601   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
602      Extract global vectors from DMDA; then duplicate for remaining
603      vectors that are the same types
604    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
605   PetscCall(DMCreateGlobalVector(da,&x));
606 
607   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
608      Set local function evaluation routine
609   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
610   switch (MMS) {
611   case 0: user.mms_solution = ZeroBCSolution; user.mms_forcing = NULL; break;
612   case 1: user.mms_solution = MMSSolution1; user.mms_forcing = MMSForcing1; break;
613   case 2: user.mms_solution = MMSSolution2; user.mms_forcing = MMSForcing2; break;
614   case 3: user.mms_solution = MMSSolution3; user.mms_forcing = MMSForcing3; break;
615   case 4: user.mms_solution = MMSSolution4; user.mms_forcing = MMSForcing4; break;
616   default: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Unknown MMS type %" PetscInt_FMT,MMS);
617   }
618 
619   if (useKokkos) {
620     PetscCheck(MMS == 1,PETSC_COMM_WORLD,PETSC_ERR_USER,"FormFunctionLocalVec_Kokkos only works with MMS 1");
621     PetscCall(DMDASNESSetFunctionLocalVec(da,INSERT_VALUES,(DMDASNESFunctionVec)FormFunctionLocalVec,&user));
622   } else {
623     PetscCall(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(DMDASNESFunction)FormFunctionLocal,&user));
624   }
625 
626   PetscCall(PetscOptionsGetBool(NULL,NULL,"-fd",&flg,NULL));
627   if (!flg) {
628     if (useKokkos) PetscCall(DMDASNESSetJacobianLocalVec(da,(DMDASNESJacobianVec)FormJacobianLocalVec,&user));
629     else PetscCall(DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,&user));
630   }
631 
632   PetscCall(PetscOptionsGetBool(NULL,NULL,"-obj",&flg,NULL));
633   if (flg) {
634     if (useKokkos) PetscCall(DMDASNESSetObjectiveLocalVec(da,(DMDASNESObjectiveVec)FormObjectiveLocalVec,&user));
635     else PetscCall(DMDASNESSetObjectiveLocal(da,(DMDASNESObjective)FormObjectiveLocal,&user));
636   }
637 
638   if (PetscDefined(HAVE_MATLAB_ENGINE)) {
639     PetscBool matlab_function = PETSC_FALSE;
640     PetscCall(PetscOptionsGetBool(NULL,NULL,"-matlab_function",&matlab_function,0));
641     if (matlab_function) {
642       PetscCall(VecDuplicate(x,&r));
643       PetscCall(SNESSetFunction(snes,r,FormFunctionMatlab,&user));
644     }
645   }
646 
647   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
648      Customize nonlinear solver; set runtime options
649    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
650   PetscCall(SNESSetFromOptions(snes));
651 
652   PetscCall(FormInitialGuess(da,&user,x));
653 
654   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
655      Solve nonlinear system
656      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
657   PetscCall(SNESSolve(snes,NULL,x));
658   PetscCall(SNESGetIterationNumber(snes,&its));
659 
660   PetscCall(SNESGetLinearSolveIterations(snes,&slits));
661   PetscCall(SNESGetKSP(snes,&ksp));
662   PetscCall(KSPGetTotalIterations(ksp,&lits));
663   PetscCheck(lits == slits,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Number of total linear iterations reported by SNES %" PetscInt_FMT " does not match reported by KSP %" PetscInt_FMT,slits,lits);
664   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
665    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
666 
667   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
668      If using MMS, check the l_2 error
669    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
670   if (setMMS) {
671     Vec       e;
672     PetscReal errorl2, errorinf;
673     PetscInt  N;
674 
675     PetscCall(VecDuplicate(x, &e));
676     PetscCall(PetscObjectViewFromOptions((PetscObject) x, NULL, "-sol_view"));
677     PetscCall(FormExactSolution(da, &user, e));
678     PetscCall(PetscObjectViewFromOptions((PetscObject) e, NULL, "-exact_view"));
679     PetscCall(VecAXPY(e, -1.0, x));
680     PetscCall(PetscObjectViewFromOptions((PetscObject) e, NULL, "-error_view"));
681     PetscCall(VecNorm(e, NORM_2, &errorl2));
682     PetscCall(VecNorm(e, NORM_INFINITY, &errorinf));
683     PetscCall(VecGetSize(e, &N));
684     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "N: %" PetscInt_FMT " error L2 %g inf %g\n", N, (double)(errorl2/PetscSqrtReal((PetscReal)N)), (double) errorinf));
685     PetscCall(VecDestroy(&e));
686     PetscCall(PetscLogEventSetDof(SNES_Solve, 0, N));
687     PetscCall(PetscLogEventSetError(SNES_Solve, 0, errorl2/PetscSqrtReal(N)));
688   }
689 
690   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
691      Free work space.  All PETSc objects should be destroyed when they
692      are no longer needed.
693    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
694   PetscCall(VecDestroy(&r));
695   PetscCall(VecDestroy(&x));
696   PetscCall(SNESDestroy(&snes));
697   PetscCall(DMDestroy(&da));
698   PetscCall(PetscFinalize());
699   return 0;
700 }
701 
702 /*TEST
703   build:
704     requires: kokkos_kernels
705     depends: ex55k.kokkos.cxx
706 
707   testset:
708     output_file: output/ex55_asm_0.out
709     requires: !single
710     args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu
711     filter: grep -v "type"
712 
713     test:
714       suffix: asm_0
715     test:
716       suffix: asm_0_kok
717       args: -use_kokkos 1 -dm_mat_type aijkokkos -dm_vec_type kokkos
718 
719 TEST*/
720