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