xref: /petsc/src/snes/tests/ex20.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1e77315bcSPatrick Sanan static char help[] = "Nonlinear Radiative Transport PDE with multigrid in 3d.\n\
2e77315bcSPatrick Sanan Uses 3-dimensional distributed arrays.\n\
3e77315bcSPatrick Sanan A 3-dim simplified Radiative Transport test problem is used, with analytic Jacobian. \n\
4e77315bcSPatrick Sanan \n\
5e77315bcSPatrick Sanan   Solves the linear systems via multilevel methods \n\
6e77315bcSPatrick Sanan \n\
7e77315bcSPatrick Sanan The command line\n\
8e77315bcSPatrick Sanan options are:\n\
9e77315bcSPatrick Sanan   -tleft <tl>, where <tl> indicates the left Diriclet BC \n\
10e77315bcSPatrick Sanan   -tright <tr>, where <tr> indicates the right Diriclet BC \n\
11e77315bcSPatrick Sanan   -beta <beta>, where <beta> indicates the exponent in T \n\n";
12e77315bcSPatrick Sanan 
13e77315bcSPatrick Sanan /*
14e77315bcSPatrick Sanan 
15e77315bcSPatrick Sanan     This example models the partial differential equation
16e77315bcSPatrick Sanan 
17e77315bcSPatrick Sanan          - Div(alpha* T^beta (GRAD T)) = 0.
18e77315bcSPatrick Sanan 
19e77315bcSPatrick Sanan     where beta = 2.5 and alpha = 1.0
20e77315bcSPatrick Sanan 
21e77315bcSPatrick Sanan     BC: T_left = 1.0, T_right = 0.1, dT/dn_top = dTdn_bottom = dT/dn_up = dT/dn_down = 0.
22e77315bcSPatrick Sanan 
23e77315bcSPatrick Sanan     in the unit square, which is uniformly discretized in each of x and
24e77315bcSPatrick Sanan     y in this simple encoding.  The degrees of freedom are cell centered.
25e77315bcSPatrick Sanan 
26e77315bcSPatrick Sanan     A finite volume approximation with the usual 7-point stencil
27e77315bcSPatrick Sanan     is used to discretize the boundary value problem to obtain a
28e77315bcSPatrick Sanan     nonlinear system of equations.
29e77315bcSPatrick Sanan 
30e77315bcSPatrick Sanan     This code was contributed by Nickolas Jovanovic based on ex18.c
31e77315bcSPatrick Sanan 
32e77315bcSPatrick Sanan */
33e77315bcSPatrick Sanan 
34e77315bcSPatrick Sanan #include <petscsnes.h>
35e77315bcSPatrick Sanan #include <petscdm.h>
36e77315bcSPatrick Sanan #include <petscdmda.h>
37e77315bcSPatrick Sanan 
38e77315bcSPatrick Sanan /* User-defined application context */
39e77315bcSPatrick Sanan 
40e77315bcSPatrick Sanan typedef struct {
41e77315bcSPatrick Sanan   PetscReal tleft, tright;   /* Dirichlet boundary conditions */
42e77315bcSPatrick Sanan   PetscReal beta, bm1, coef; /* nonlinear diffusivity parameterizations */
43379ef8d2SAlexander   PetscBool converged;       /* For user convergence test, whether to mark as converged */
44e77315bcSPatrick Sanan } AppCtx;
45e77315bcSPatrick Sanan 
46e77315bcSPatrick Sanan #define POWFLOP 5 /* assume a pow() takes five flops */
47e77315bcSPatrick Sanan 
48e77315bcSPatrick Sanan extern PetscErrorCode FormInitialGuess(SNES, Vec, void *);
49e77315bcSPatrick Sanan extern PetscErrorCode FormFunction(SNES, Vec, Vec, void *);
50e77315bcSPatrick Sanan extern PetscErrorCode FormJacobian(SNES, Vec, Mat, Mat, void *);
51379ef8d2SAlexander extern PetscErrorCode TestConvergence(SNES, PetscInt, PetscReal, PetscReal, PetscReal, SNESConvergedReason *, void *);
52e77315bcSPatrick Sanan 
main(int argc,char ** argv)53d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
54d71ae5a4SJacob Faibussowitsch {
55e77315bcSPatrick Sanan   SNES      snes;
56e77315bcSPatrick Sanan   AppCtx    user;
57e77315bcSPatrick Sanan   PetscInt  its, lits;
58379ef8d2SAlexander   PetscReal litspit = 0; /* avoid uninitialized warning */
59e77315bcSPatrick Sanan   DM        da;
60379ef8d2SAlexander   PetscBool use_convergence_test = PETSC_FALSE;
61e77315bcSPatrick Sanan 
62327415f7SBarry Smith   PetscFunctionBeginUser;
639566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
64e77315bcSPatrick Sanan   /* set problem parameters */
65e77315bcSPatrick Sanan   user.tleft  = 1.0;
66e77315bcSPatrick Sanan   user.tright = 0.1;
67e77315bcSPatrick Sanan   user.beta   = 2.5;
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-tleft", &user.tleft, NULL));
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-tright", &user.tright, NULL));
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-beta", &user.beta, NULL));
71379ef8d2SAlexander   PetscCall(PetscOptionsGetBool(NULL, NULL, "-use_convergence_test", &use_convergence_test, NULL));
72379ef8d2SAlexander   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mark_converged", &user.converged, NULL));
73e77315bcSPatrick Sanan   user.bm1  = user.beta - 1.0;
74e77315bcSPatrick Sanan   user.coef = user.beta / 2.0;
75e77315bcSPatrick Sanan 
76e77315bcSPatrick Sanan   /*
77e77315bcSPatrick Sanan       Set the DMDA (grid structure) for the grids.
78e77315bcSPatrick Sanan   */
799566063dSJacob Faibussowitsch   PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 5, 5, 5, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 1, 1, 0, 0, 0, &da));
809566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
819566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
829566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, &user));
83e77315bcSPatrick Sanan 
84e77315bcSPatrick Sanan   /*
85e77315bcSPatrick Sanan      Create the nonlinear solver
86e77315bcSPatrick Sanan   */
879566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
889566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes, da));
899566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, FormFunction, &user));
909566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes, NULL, NULL, FormJacobian, &user));
91379ef8d2SAlexander   if (use_convergence_test) PetscCall(SNESSetConvergenceTest(snes, TestConvergence, &user, NULL));
929566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
939566063dSJacob Faibussowitsch   PetscCall(SNESSetComputeInitialGuess(snes, FormInitialGuess, NULL));
94e77315bcSPatrick Sanan 
959566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes, NULL, NULL));
969566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes, &its));
979566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(snes, &lits));
98379ef8d2SAlexander   if (its) litspit = ((PetscReal)lits) / ((PetscReal)its);
9963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of SNES iterations = %" PetscInt_FMT "\n", its));
10063a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of Linear iterations = %" PetscInt_FMT "\n", lits));
101379ef8d2SAlexander   if (its) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Average Linear its / SNES = %e\n", (double)litspit));
102e77315bcSPatrick Sanan 
1039566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1049566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1059566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
106b122ec5aSJacob Faibussowitsch   return 0;
107e77315bcSPatrick Sanan }
108e77315bcSPatrick Sanan /* --------------------  Form initial approximation ----------------- */
FormInitialGuess(SNES snes,Vec X,PetscCtx ctx)109*2a8381b2SBarry Smith PetscErrorCode FormInitialGuess(SNES snes, Vec X, PetscCtx ctx)
110d71ae5a4SJacob Faibussowitsch {
111e77315bcSPatrick Sanan   AppCtx        *user;
112e77315bcSPatrick Sanan   PetscInt       i, j, k, xs, ys, xm, ym, zs, zm;
113e77315bcSPatrick Sanan   PetscScalar ***x;
114e77315bcSPatrick Sanan   DM             da;
115e77315bcSPatrick Sanan 
116e77315bcSPatrick Sanan   PetscFunctionBeginUser;
1179566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
1189566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da, &user));
1199566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1209566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
121e77315bcSPatrick Sanan 
122e77315bcSPatrick Sanan   /* Compute initial guess */
123e77315bcSPatrick Sanan   for (k = zs; k < zs + zm; k++) {
124e77315bcSPatrick Sanan     for (j = ys; j < ys + ym; j++) {
125ad540459SPierre Jolivet       for (i = xs; i < xs + xm; i++) x[k][j][i] = user->tleft;
126e77315bcSPatrick Sanan     }
127e77315bcSPatrick Sanan   }
1289566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
1293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
130e77315bcSPatrick Sanan }
131e77315bcSPatrick Sanan /* --------------------  Evaluate Function F(x) --------------------- */
FormFunction(SNES snes,Vec X,Vec F,void * ptr)132d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *ptr)
133d71ae5a4SJacob Faibussowitsch {
134e77315bcSPatrick Sanan   AppCtx        *user = (AppCtx *)ptr;
135e77315bcSPatrick Sanan   PetscInt       i, j, k, mx, my, mz, xs, ys, zs, xm, ym, zm;
136e77315bcSPatrick Sanan   PetscScalar    zero = 0.0, one = 1.0;
137e77315bcSPatrick Sanan   PetscScalar    hx, hy, hz, hxhydhz, hyhzdhx, hzhxdhy;
138e77315bcSPatrick Sanan   PetscScalar    t0, tn, ts, te, tw, an, as, ae, aw, dn, ds, de, dw, fn = 0.0, fs = 0.0, fe = 0.0, fw = 0.0;
139e77315bcSPatrick Sanan   PetscScalar    tleft, tright, beta, td, ad, dd, fd = 0.0, tu, au, du = 0.0, fu = 0.0;
140e77315bcSPatrick Sanan   PetscScalar ***x, ***f;
141e77315bcSPatrick Sanan   Vec            localX;
142e77315bcSPatrick Sanan   DM             da;
143e77315bcSPatrick Sanan 
144e77315bcSPatrick Sanan   PetscFunctionBeginUser;
1459566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
1469566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
1479566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
1489371c9d4SSatish Balay   hx      = one / (PetscReal)(mx - 1);
1499371c9d4SSatish Balay   hy      = one / (PetscReal)(my - 1);
1509371c9d4SSatish Balay   hz      = one / (PetscReal)(mz - 1);
1519371c9d4SSatish Balay   hxhydhz = hx * hy / hz;
1529371c9d4SSatish Balay   hyhzdhx = hy * hz / hx;
1539371c9d4SSatish Balay   hzhxdhy = hz * hx / hy;
1549371c9d4SSatish Balay   tleft   = user->tleft;
1559371c9d4SSatish Balay   tright  = user->tright;
156e77315bcSPatrick Sanan   beta    = user->beta;
157e77315bcSPatrick Sanan 
158e77315bcSPatrick Sanan   /* Get ghost points */
1599566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
1609566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
1619566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1629566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, localX, &x));
1639566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
164e77315bcSPatrick Sanan 
165e77315bcSPatrick Sanan   /* Evaluate function */
166e77315bcSPatrick Sanan   for (k = zs; k < zs + zm; k++) {
167e77315bcSPatrick Sanan     for (j = ys; j < ys + ym; j++) {
168e77315bcSPatrick Sanan       for (i = xs; i < xs + xm; i++) {
169e77315bcSPatrick Sanan         t0 = x[k][j][i];
170e77315bcSPatrick Sanan 
171e77315bcSPatrick Sanan         if (i > 0 && i < mx - 1 && j > 0 && j < my - 1 && k > 0 && k < mz - 1) {
172e77315bcSPatrick Sanan           /* general interior volume */
173e77315bcSPatrick Sanan 
174e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
175e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
176e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
177e77315bcSPatrick Sanan           fw = dw * (t0 - tw);
178e77315bcSPatrick Sanan 
179e77315bcSPatrick Sanan           te = x[k][j][i + 1];
180e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
181e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
182e77315bcSPatrick Sanan           fe = de * (te - t0);
183e77315bcSPatrick Sanan 
184e77315bcSPatrick Sanan           ts = x[k][j - 1][i];
185e77315bcSPatrick Sanan           as = 0.5 * (t0 + ts);
186e77315bcSPatrick Sanan           ds = PetscPowScalar(as, beta);
187e77315bcSPatrick Sanan           fs = ds * (t0 - ts);
188e77315bcSPatrick Sanan 
189e77315bcSPatrick Sanan           tn = x[k][j + 1][i];
190e77315bcSPatrick Sanan           an = 0.5 * (t0 + tn);
191e77315bcSPatrick Sanan           dn = PetscPowScalar(an, beta);
192e77315bcSPatrick Sanan           fn = dn * (tn - t0);
193e77315bcSPatrick Sanan 
194e77315bcSPatrick Sanan           td = x[k - 1][j][i];
195e77315bcSPatrick Sanan           ad = 0.5 * (t0 + td);
196e77315bcSPatrick Sanan           dd = PetscPowScalar(ad, beta);
197e77315bcSPatrick Sanan           fd = dd * (t0 - td);
198e77315bcSPatrick Sanan 
199e77315bcSPatrick Sanan           tu = x[k + 1][j][i];
200e77315bcSPatrick Sanan           au = 0.5 * (t0 + tu);
201e77315bcSPatrick Sanan           du = PetscPowScalar(au, beta);
202e77315bcSPatrick Sanan           fu = du * (tu - t0);
203e77315bcSPatrick Sanan 
204e77315bcSPatrick Sanan         } else if (i == 0) {
205e77315bcSPatrick Sanan           /* left-hand (west) boundary */
206e77315bcSPatrick Sanan           tw = tleft;
207e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
208e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
209e77315bcSPatrick Sanan           fw = dw * (t0 - tw);
210e77315bcSPatrick Sanan 
211e77315bcSPatrick Sanan           te = x[k][j][i + 1];
212e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
213e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
214e77315bcSPatrick Sanan           fe = de * (te - t0);
215e77315bcSPatrick Sanan 
216e77315bcSPatrick Sanan           if (j > 0) {
217e77315bcSPatrick Sanan             ts = x[k][j - 1][i];
218e77315bcSPatrick Sanan             as = 0.5 * (t0 + ts);
219e77315bcSPatrick Sanan             ds = PetscPowScalar(as, beta);
220e77315bcSPatrick Sanan             fs = ds * (t0 - ts);
221e77315bcSPatrick Sanan           } else {
222e77315bcSPatrick Sanan             fs = zero;
223e77315bcSPatrick Sanan           }
224e77315bcSPatrick Sanan 
225e77315bcSPatrick Sanan           if (j < my - 1) {
226e77315bcSPatrick Sanan             tn = x[k][j + 1][i];
227e77315bcSPatrick Sanan             an = 0.5 * (t0 + tn);
228e77315bcSPatrick Sanan             dn = PetscPowScalar(an, beta);
229e77315bcSPatrick Sanan             fn = dn * (tn - t0);
230e77315bcSPatrick Sanan           } else {
231e77315bcSPatrick Sanan             fn = zero;
232e77315bcSPatrick Sanan           }
233e77315bcSPatrick Sanan 
234e77315bcSPatrick Sanan           if (k > 0) {
235e77315bcSPatrick Sanan             td = x[k - 1][j][i];
236e77315bcSPatrick Sanan             ad = 0.5 * (t0 + td);
237e77315bcSPatrick Sanan             dd = PetscPowScalar(ad, beta);
238e77315bcSPatrick Sanan             fd = dd * (t0 - td);
239e77315bcSPatrick Sanan           } else {
240e77315bcSPatrick Sanan             fd = zero;
241e77315bcSPatrick Sanan           }
242e77315bcSPatrick Sanan 
243e77315bcSPatrick Sanan           if (k < mz - 1) {
244e77315bcSPatrick Sanan             tu = x[k + 1][j][i];
245e77315bcSPatrick Sanan             au = 0.5 * (t0 + tu);
246e77315bcSPatrick Sanan             du = PetscPowScalar(au, beta);
247e77315bcSPatrick Sanan             fu = du * (tu - t0);
248e77315bcSPatrick Sanan           } else {
249e77315bcSPatrick Sanan             fu = zero;
250e77315bcSPatrick Sanan           }
251e77315bcSPatrick Sanan 
252e77315bcSPatrick Sanan         } else if (i == mx - 1) {
253e77315bcSPatrick Sanan           /* right-hand (east) boundary */
254e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
255e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
256e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
257e77315bcSPatrick Sanan           fw = dw * (t0 - tw);
258e77315bcSPatrick Sanan 
259e77315bcSPatrick Sanan           te = tright;
260e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
261e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
262e77315bcSPatrick Sanan           fe = de * (te - t0);
263e77315bcSPatrick Sanan 
264e77315bcSPatrick Sanan           if (j > 0) {
265e77315bcSPatrick Sanan             ts = x[k][j - 1][i];
266e77315bcSPatrick Sanan             as = 0.5 * (t0 + ts);
267e77315bcSPatrick Sanan             ds = PetscPowScalar(as, beta);
268e77315bcSPatrick Sanan             fs = ds * (t0 - ts);
269e77315bcSPatrick Sanan           } else {
270e77315bcSPatrick Sanan             fs = zero;
271e77315bcSPatrick Sanan           }
272e77315bcSPatrick Sanan 
273e77315bcSPatrick Sanan           if (j < my - 1) {
274e77315bcSPatrick Sanan             tn = x[k][j + 1][i];
275e77315bcSPatrick Sanan             an = 0.5 * (t0 + tn);
276e77315bcSPatrick Sanan             dn = PetscPowScalar(an, beta);
277e77315bcSPatrick Sanan             fn = dn * (tn - t0);
278e77315bcSPatrick Sanan           } else {
279e77315bcSPatrick Sanan             fn = zero;
280e77315bcSPatrick Sanan           }
281e77315bcSPatrick Sanan 
282e77315bcSPatrick Sanan           if (k > 0) {
283e77315bcSPatrick Sanan             td = x[k - 1][j][i];
284e77315bcSPatrick Sanan             ad = 0.5 * (t0 + td);
285e77315bcSPatrick Sanan             dd = PetscPowScalar(ad, beta);
286e77315bcSPatrick Sanan             fd = dd * (t0 - td);
287e77315bcSPatrick Sanan           } else {
288e77315bcSPatrick Sanan             fd = zero;
289e77315bcSPatrick Sanan           }
290e77315bcSPatrick Sanan 
291e77315bcSPatrick Sanan           if (k < mz - 1) {
292e77315bcSPatrick Sanan             tu = x[k + 1][j][i];
293e77315bcSPatrick Sanan             au = 0.5 * (t0 + tu);
294e77315bcSPatrick Sanan             du = PetscPowScalar(au, beta);
295e77315bcSPatrick Sanan             fu = du * (tu - t0);
296e77315bcSPatrick Sanan           } else {
297e77315bcSPatrick Sanan             fu = zero;
298e77315bcSPatrick Sanan           }
299e77315bcSPatrick Sanan 
300e77315bcSPatrick Sanan         } else if (j == 0) {
301e77315bcSPatrick Sanan           /* bottom (south) boundary, and i <> 0 or mx-1 */
302e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
303e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
304e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
305e77315bcSPatrick Sanan           fw = dw * (t0 - tw);
306e77315bcSPatrick Sanan 
307e77315bcSPatrick Sanan           te = x[k][j][i + 1];
308e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
309e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
310e77315bcSPatrick Sanan           fe = de * (te - t0);
311e77315bcSPatrick Sanan 
312e77315bcSPatrick Sanan           fs = zero;
313e77315bcSPatrick Sanan 
314e77315bcSPatrick Sanan           tn = x[k][j + 1][i];
315e77315bcSPatrick Sanan           an = 0.5 * (t0 + tn);
316e77315bcSPatrick Sanan           dn = PetscPowScalar(an, beta);
317e77315bcSPatrick Sanan           fn = dn * (tn - t0);
318e77315bcSPatrick Sanan 
319e77315bcSPatrick Sanan           if (k > 0) {
320e77315bcSPatrick Sanan             td = x[k - 1][j][i];
321e77315bcSPatrick Sanan             ad = 0.5 * (t0 + td);
322e77315bcSPatrick Sanan             dd = PetscPowScalar(ad, beta);
323e77315bcSPatrick Sanan             fd = dd * (t0 - td);
324e77315bcSPatrick Sanan           } else {
325e77315bcSPatrick Sanan             fd = zero;
326e77315bcSPatrick Sanan           }
327e77315bcSPatrick Sanan 
328e77315bcSPatrick Sanan           if (k < mz - 1) {
329e77315bcSPatrick Sanan             tu = x[k + 1][j][i];
330e77315bcSPatrick Sanan             au = 0.5 * (t0 + tu);
331e77315bcSPatrick Sanan             du = PetscPowScalar(au, beta);
332e77315bcSPatrick Sanan             fu = du * (tu - t0);
333e77315bcSPatrick Sanan           } else {
334e77315bcSPatrick Sanan             fu = zero;
335e77315bcSPatrick Sanan           }
336e77315bcSPatrick Sanan 
337e77315bcSPatrick Sanan         } else if (j == my - 1) {
338e77315bcSPatrick Sanan           /* top (north) boundary, and i <> 0 or mx-1 */
339e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
340e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
341e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
342e77315bcSPatrick Sanan           fw = dw * (t0 - tw);
343e77315bcSPatrick Sanan 
344e77315bcSPatrick Sanan           te = x[k][j][i + 1];
345e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
346e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
347e77315bcSPatrick Sanan           fe = de * (te - t0);
348e77315bcSPatrick Sanan 
349e77315bcSPatrick Sanan           ts = x[k][j - 1][i];
350e77315bcSPatrick Sanan           as = 0.5 * (t0 + ts);
351e77315bcSPatrick Sanan           ds = PetscPowScalar(as, beta);
352e77315bcSPatrick Sanan           fs = ds * (t0 - ts);
353e77315bcSPatrick Sanan 
354e77315bcSPatrick Sanan           fn = zero;
355e77315bcSPatrick Sanan 
356e77315bcSPatrick Sanan           if (k > 0) {
357e77315bcSPatrick Sanan             td = x[k - 1][j][i];
358e77315bcSPatrick Sanan             ad = 0.5 * (t0 + td);
359e77315bcSPatrick Sanan             dd = PetscPowScalar(ad, beta);
360e77315bcSPatrick Sanan             fd = dd * (t0 - td);
361e77315bcSPatrick Sanan           } else {
362e77315bcSPatrick Sanan             fd = zero;
363e77315bcSPatrick Sanan           }
364e77315bcSPatrick Sanan 
365e77315bcSPatrick Sanan           if (k < mz - 1) {
366e77315bcSPatrick Sanan             tu = x[k + 1][j][i];
367e77315bcSPatrick Sanan             au = 0.5 * (t0 + tu);
368e77315bcSPatrick Sanan             du = PetscPowScalar(au, beta);
369e77315bcSPatrick Sanan             fu = du * (tu - t0);
370e77315bcSPatrick Sanan           } else {
371e77315bcSPatrick Sanan             fu = zero;
372e77315bcSPatrick Sanan           }
373e77315bcSPatrick Sanan 
374e77315bcSPatrick Sanan         } else if (k == 0) {
375e77315bcSPatrick Sanan           /* down boundary (interior only) */
376e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
377e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
378e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
379e77315bcSPatrick Sanan           fw = dw * (t0 - tw);
380e77315bcSPatrick Sanan 
381e77315bcSPatrick Sanan           te = x[k][j][i + 1];
382e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
383e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
384e77315bcSPatrick Sanan           fe = de * (te - t0);
385e77315bcSPatrick Sanan 
386e77315bcSPatrick Sanan           ts = x[k][j - 1][i];
387e77315bcSPatrick Sanan           as = 0.5 * (t0 + ts);
388e77315bcSPatrick Sanan           ds = PetscPowScalar(as, beta);
389e77315bcSPatrick Sanan           fs = ds * (t0 - ts);
390e77315bcSPatrick Sanan 
391e77315bcSPatrick Sanan           tn = x[k][j + 1][i];
392e77315bcSPatrick Sanan           an = 0.5 * (t0 + tn);
393e77315bcSPatrick Sanan           dn = PetscPowScalar(an, beta);
394e77315bcSPatrick Sanan           fn = dn * (tn - t0);
395e77315bcSPatrick Sanan 
396e77315bcSPatrick Sanan           fd = zero;
397e77315bcSPatrick Sanan 
398e77315bcSPatrick Sanan           tu = x[k + 1][j][i];
399e77315bcSPatrick Sanan           au = 0.5 * (t0 + tu);
400e77315bcSPatrick Sanan           du = PetscPowScalar(au, beta);
401e77315bcSPatrick Sanan           fu = du * (tu - t0);
402e77315bcSPatrick Sanan 
403e77315bcSPatrick Sanan         } else if (k == mz - 1) {
404e77315bcSPatrick Sanan           /* up boundary (interior only) */
405e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
406e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
407e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
408e77315bcSPatrick Sanan           fw = dw * (t0 - tw);
409e77315bcSPatrick Sanan 
410e77315bcSPatrick Sanan           te = x[k][j][i + 1];
411e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
412e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
413e77315bcSPatrick Sanan           fe = de * (te - t0);
414e77315bcSPatrick Sanan 
415e77315bcSPatrick Sanan           ts = x[k][j - 1][i];
416e77315bcSPatrick Sanan           as = 0.5 * (t0 + ts);
417e77315bcSPatrick Sanan           ds = PetscPowScalar(as, beta);
418e77315bcSPatrick Sanan           fs = ds * (t0 - ts);
419e77315bcSPatrick Sanan 
420e77315bcSPatrick Sanan           tn = x[k][j + 1][i];
421e77315bcSPatrick Sanan           an = 0.5 * (t0 + tn);
422e77315bcSPatrick Sanan           dn = PetscPowScalar(an, beta);
423e77315bcSPatrick Sanan           fn = dn * (tn - t0);
424e77315bcSPatrick Sanan 
425e77315bcSPatrick Sanan           td = x[k - 1][j][i];
426e77315bcSPatrick Sanan           ad = 0.5 * (t0 + td);
427e77315bcSPatrick Sanan           dd = PetscPowScalar(ad, beta);
428e77315bcSPatrick Sanan           fd = dd * (t0 - td);
429e77315bcSPatrick Sanan 
430e77315bcSPatrick Sanan           fu = zero;
431e77315bcSPatrick Sanan         }
432e77315bcSPatrick Sanan 
433e77315bcSPatrick Sanan         f[k][j][i] = -hyhzdhx * (fe - fw) - hzhxdhy * (fn - fs) - hxhydhz * (fu - fd);
434e77315bcSPatrick Sanan       }
435e77315bcSPatrick Sanan     }
436e77315bcSPatrick Sanan   }
4379566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, localX, &x));
4389566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
4399566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
4409566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops((22.0 + 4.0 * POWFLOP) * ym * xm));
4413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
442e77315bcSPatrick Sanan }
443e77315bcSPatrick Sanan /* --------------------  Evaluate Jacobian F(x) --------------------- */
FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void * ptr)444d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat jac, void *ptr)
445d71ae5a4SJacob Faibussowitsch {
446e77315bcSPatrick Sanan   AppCtx        *user = (AppCtx *)ptr;
447e77315bcSPatrick Sanan   PetscInt       i, j, k, mx, my, mz, xs, ys, zs, xm, ym, zm;
448e77315bcSPatrick Sanan   PetscScalar    one = 1.0;
449e77315bcSPatrick Sanan   PetscScalar    hx, hy, hz, hxhydhz, hyhzdhx, hzhxdhy;
450e77315bcSPatrick Sanan   PetscScalar    t0, tn, ts, te, tw, an, as, ae, aw, dn, ds, de, dw;
451e77315bcSPatrick Sanan   PetscScalar    tleft, tright, beta, td, ad, dd, tu, au, du, v[7], bm1, coef;
452e77315bcSPatrick Sanan   PetscScalar ***x, bn, bs, be, bw, bu, bd, gn, gs, ge, gw, gu, gd;
453e77315bcSPatrick Sanan   Vec            localX;
454e77315bcSPatrick Sanan   MatStencil     c[7], row;
455e77315bcSPatrick Sanan   DM             da;
456e77315bcSPatrick Sanan 
457e77315bcSPatrick Sanan   PetscFunctionBeginUser;
4589566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &da));
4599566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
4609566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, NULL, &mx, &my, &mz, 0, 0, 0, 0, 0, 0, 0, 0, 0));
4619371c9d4SSatish Balay   hx      = one / (PetscReal)(mx - 1);
4629371c9d4SSatish Balay   hy      = one / (PetscReal)(my - 1);
4639371c9d4SSatish Balay   hz      = one / (PetscReal)(mz - 1);
4649371c9d4SSatish Balay   hxhydhz = hx * hy / hz;
4659371c9d4SSatish Balay   hyhzdhx = hy * hz / hx;
4669371c9d4SSatish Balay   hzhxdhy = hz * hx / hy;
4679371c9d4SSatish Balay   tleft   = user->tleft;
4689371c9d4SSatish Balay   tright  = user->tright;
4699371c9d4SSatish Balay   beta    = user->beta;
4709371c9d4SSatish Balay   bm1     = user->bm1;
4719371c9d4SSatish Balay   coef    = user->coef;
472e77315bcSPatrick Sanan 
473e77315bcSPatrick Sanan   /* Get ghost points */
4749566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
4759566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
4769566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
4779566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, localX, &x));
478e77315bcSPatrick Sanan 
479e77315bcSPatrick Sanan   /* Evaluate Jacobian of function */
480e77315bcSPatrick Sanan   for (k = zs; k < zs + zm; k++) {
481e77315bcSPatrick Sanan     for (j = ys; j < ys + ym; j++) {
482e77315bcSPatrick Sanan       for (i = xs; i < xs + xm; i++) {
483e77315bcSPatrick Sanan         t0    = x[k][j][i];
4849371c9d4SSatish Balay         row.k = k;
4859371c9d4SSatish Balay         row.j = j;
4869371c9d4SSatish Balay         row.i = i;
487e77315bcSPatrick Sanan         if (i > 0 && i < mx - 1 && j > 0 && j < my - 1 && k > 0 && k < mz - 1) {
488e77315bcSPatrick Sanan           /* general interior volume */
489e77315bcSPatrick Sanan 
490e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
491e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
492e77315bcSPatrick Sanan           bw = PetscPowScalar(aw, bm1);
493e77315bcSPatrick Sanan           /* dw = bw * aw */
494e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
495e77315bcSPatrick Sanan           gw = coef * bw * (t0 - tw);
496e77315bcSPatrick Sanan 
497e77315bcSPatrick Sanan           te = x[k][j][i + 1];
498e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
499e77315bcSPatrick Sanan           be = PetscPowScalar(ae, bm1);
500e77315bcSPatrick Sanan           /* de = be * ae; */
501e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
502e77315bcSPatrick Sanan           ge = coef * be * (te - t0);
503e77315bcSPatrick Sanan 
504e77315bcSPatrick Sanan           ts = x[k][j - 1][i];
505e77315bcSPatrick Sanan           as = 0.5 * (t0 + ts);
506e77315bcSPatrick Sanan           bs = PetscPowScalar(as, bm1);
507e77315bcSPatrick Sanan           /* ds = bs * as; */
508e77315bcSPatrick Sanan           ds = PetscPowScalar(as, beta);
509e77315bcSPatrick Sanan           gs = coef * bs * (t0 - ts);
510e77315bcSPatrick Sanan 
511e77315bcSPatrick Sanan           tn = x[k][j + 1][i];
512e77315bcSPatrick Sanan           an = 0.5 * (t0 + tn);
513e77315bcSPatrick Sanan           bn = PetscPowScalar(an, bm1);
514e77315bcSPatrick Sanan           /* dn = bn * an; */
515e77315bcSPatrick Sanan           dn = PetscPowScalar(an, beta);
516e77315bcSPatrick Sanan           gn = coef * bn * (tn - t0);
517e77315bcSPatrick Sanan 
518e77315bcSPatrick Sanan           td = x[k - 1][j][i];
519e77315bcSPatrick Sanan           ad = 0.5 * (t0 + td);
520e77315bcSPatrick Sanan           bd = PetscPowScalar(ad, bm1);
521e77315bcSPatrick Sanan           /* dd = bd * ad; */
522e77315bcSPatrick Sanan           dd = PetscPowScalar(ad, beta);
523e77315bcSPatrick Sanan           gd = coef * bd * (t0 - td);
524e77315bcSPatrick Sanan 
525e77315bcSPatrick Sanan           tu = x[k + 1][j][i];
526e77315bcSPatrick Sanan           au = 0.5 * (t0 + tu);
527e77315bcSPatrick Sanan           bu = PetscPowScalar(au, bm1);
528e77315bcSPatrick Sanan           /* du = bu * au; */
529e77315bcSPatrick Sanan           du = PetscPowScalar(au, beta);
530e77315bcSPatrick Sanan           gu = coef * bu * (tu - t0);
531e77315bcSPatrick Sanan 
5329371c9d4SSatish Balay           c[0].k = k - 1;
5339371c9d4SSatish Balay           c[0].j = j;
5349371c9d4SSatish Balay           c[0].i = i;
5359371c9d4SSatish Balay           v[0]   = -hxhydhz * (dd - gd);
5369371c9d4SSatish Balay           c[1].k = k;
5379371c9d4SSatish Balay           c[1].j = j - 1;
5389371c9d4SSatish Balay           c[1].i = i;
539e77315bcSPatrick Sanan           v[1]   = -hzhxdhy * (ds - gs);
5409371c9d4SSatish Balay           c[2].k = k;
5419371c9d4SSatish Balay           c[2].j = j;
5429371c9d4SSatish Balay           c[2].i = i - 1;
543e77315bcSPatrick Sanan           v[2]   = -hyhzdhx * (dw - gw);
5449371c9d4SSatish Balay           c[3].k = k;
5459371c9d4SSatish Balay           c[3].j = j;
5469371c9d4SSatish Balay           c[3].i = i;
547e77315bcSPatrick Sanan           v[3]   = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu);
5489371c9d4SSatish Balay           c[4].k = k;
5499371c9d4SSatish Balay           c[4].j = j;
5509371c9d4SSatish Balay           c[4].i = i + 1;
551e77315bcSPatrick Sanan           v[4]   = -hyhzdhx * (de + ge);
5529371c9d4SSatish Balay           c[5].k = k;
5539371c9d4SSatish Balay           c[5].j = j + 1;
5549371c9d4SSatish Balay           c[5].i = i;
555e77315bcSPatrick Sanan           v[5]   = -hzhxdhy * (dn + gn);
5569371c9d4SSatish Balay           c[6].k = k + 1;
5579371c9d4SSatish Balay           c[6].j = j;
5589371c9d4SSatish Balay           c[6].i = i;
559e77315bcSPatrick Sanan           v[6]   = -hxhydhz * (du + gu);
5609566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(jac, 1, &row, 7, c, v, INSERT_VALUES));
561e77315bcSPatrick Sanan 
562e77315bcSPatrick Sanan         } else if (i == 0) {
563e77315bcSPatrick Sanan           /* left-hand plane boundary */
564e77315bcSPatrick Sanan           tw = tleft;
565e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
566e77315bcSPatrick Sanan           bw = PetscPowScalar(aw, bm1);
567e77315bcSPatrick Sanan           /* dw = bw * aw */
568e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
569e77315bcSPatrick Sanan           gw = coef * bw * (t0 - tw);
570e77315bcSPatrick Sanan 
571e77315bcSPatrick Sanan           te = x[k][j][i + 1];
572e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
573e77315bcSPatrick Sanan           be = PetscPowScalar(ae, bm1);
574e77315bcSPatrick Sanan           /* de = be * ae; */
575e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
576e77315bcSPatrick Sanan           ge = coef * be * (te - t0);
577e77315bcSPatrick Sanan 
578e77315bcSPatrick Sanan           /* left-hand bottom edge */
579e77315bcSPatrick Sanan           if (j == 0) {
580e77315bcSPatrick Sanan             tn = x[k][j + 1][i];
581e77315bcSPatrick Sanan             an = 0.5 * (t0 + tn);
582e77315bcSPatrick Sanan             bn = PetscPowScalar(an, bm1);
583e77315bcSPatrick Sanan             /* dn = bn * an; */
584e77315bcSPatrick Sanan             dn = PetscPowScalar(an, beta);
585e77315bcSPatrick Sanan             gn = coef * bn * (tn - t0);
586e77315bcSPatrick Sanan 
587e77315bcSPatrick Sanan             /* left-hand bottom down corner */
588e77315bcSPatrick Sanan             if (k == 0) {
589e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
590e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
591e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
592e77315bcSPatrick Sanan               /* du = bu * au; */
593e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
594e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
595e77315bcSPatrick Sanan 
5969371c9d4SSatish Balay               c[0].k = k;
5979371c9d4SSatish Balay               c[0].j = j;
5989371c9d4SSatish Balay               c[0].i = i;
599e77315bcSPatrick Sanan               v[0]   = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu);
6009371c9d4SSatish Balay               c[1].k = k;
6019371c9d4SSatish Balay               c[1].j = j;
6029371c9d4SSatish Balay               c[1].i = i + 1;
603e77315bcSPatrick Sanan               v[1]   = -hyhzdhx * (de + ge);
6049371c9d4SSatish Balay               c[2].k = k;
6059371c9d4SSatish Balay               c[2].j = j + 1;
6069371c9d4SSatish Balay               c[2].i = i;
607e77315bcSPatrick Sanan               v[2]   = -hzhxdhy * (dn + gn);
6089371c9d4SSatish Balay               c[3].k = k + 1;
6099371c9d4SSatish Balay               c[3].j = j;
6109371c9d4SSatish Balay               c[3].i = i;
611e77315bcSPatrick Sanan               v[3]   = -hxhydhz * (du + gu);
6129566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES));
613e77315bcSPatrick Sanan 
614e77315bcSPatrick Sanan               /* left-hand bottom interior edge */
615e77315bcSPatrick Sanan             } else if (k < mz - 1) {
616e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
617e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
618e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
619e77315bcSPatrick Sanan               /* du = bu * au; */
620e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
621e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
622e77315bcSPatrick Sanan 
623e77315bcSPatrick Sanan               td = x[k - 1][j][i];
624e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
625e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
626e77315bcSPatrick Sanan               /* dd = bd * ad; */
627e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
628e77315bcSPatrick Sanan               gd = coef * bd * (td - t0);
629e77315bcSPatrick Sanan 
6309371c9d4SSatish Balay               c[0].k = k - 1;
6319371c9d4SSatish Balay               c[0].j = j;
6329371c9d4SSatish Balay               c[0].i = i;
633e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
6349371c9d4SSatish Balay               c[1].k = k;
6359371c9d4SSatish Balay               c[1].j = j;
6369371c9d4SSatish Balay               c[1].i = i;
637e77315bcSPatrick Sanan               v[1]   = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu);
6389371c9d4SSatish Balay               c[2].k = k;
6399371c9d4SSatish Balay               c[2].j = j;
6409371c9d4SSatish Balay               c[2].i = i + 1;
641e77315bcSPatrick Sanan               v[2]   = -hyhzdhx * (de + ge);
6429371c9d4SSatish Balay               c[3].k = k;
6439371c9d4SSatish Balay               c[3].j = j + 1;
6449371c9d4SSatish Balay               c[3].i = i;
645e77315bcSPatrick Sanan               v[3]   = -hzhxdhy * (dn + gn);
6469371c9d4SSatish Balay               c[4].k = k + 1;
6479371c9d4SSatish Balay               c[4].j = j;
6489371c9d4SSatish Balay               c[4].i = i;
649e77315bcSPatrick Sanan               v[4]   = -hxhydhz * (du + gu);
6509566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
651e77315bcSPatrick Sanan 
652e77315bcSPatrick Sanan               /* left-hand bottom up corner */
653e77315bcSPatrick Sanan             } else {
654e77315bcSPatrick Sanan               td = x[k - 1][j][i];
655e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
656e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
657e77315bcSPatrick Sanan               /* dd = bd * ad; */
658e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
659e77315bcSPatrick Sanan               gd = coef * bd * (td - t0);
660e77315bcSPatrick Sanan 
6619371c9d4SSatish Balay               c[0].k = k - 1;
6629371c9d4SSatish Balay               c[0].j = j;
6639371c9d4SSatish Balay               c[0].i = i;
664e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
6659371c9d4SSatish Balay               c[1].k = k;
6669371c9d4SSatish Balay               c[1].j = j;
6679371c9d4SSatish Balay               c[1].i = i;
668e77315bcSPatrick Sanan               v[1]   = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd);
6699371c9d4SSatish Balay               c[2].k = k;
6709371c9d4SSatish Balay               c[2].j = j;
6719371c9d4SSatish Balay               c[2].i = i + 1;
672e77315bcSPatrick Sanan               v[2]   = -hyhzdhx * (de + ge);
6739371c9d4SSatish Balay               c[3].k = k;
6749371c9d4SSatish Balay               c[3].j = j + 1;
6759371c9d4SSatish Balay               c[3].i = i;
676e77315bcSPatrick Sanan               v[3]   = -hzhxdhy * (dn + gn);
6779566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES));
678e77315bcSPatrick Sanan             }
679e77315bcSPatrick Sanan 
680e77315bcSPatrick Sanan             /* left-hand top edge */
681e77315bcSPatrick Sanan           } else if (j == my - 1) {
682e77315bcSPatrick Sanan             ts = x[k][j - 1][i];
683e77315bcSPatrick Sanan             as = 0.5 * (t0 + ts);
684e77315bcSPatrick Sanan             bs = PetscPowScalar(as, bm1);
685e77315bcSPatrick Sanan             /* ds = bs * as; */
686e77315bcSPatrick Sanan             ds = PetscPowScalar(as, beta);
687e77315bcSPatrick Sanan             gs = coef * bs * (ts - t0);
688e77315bcSPatrick Sanan 
689e77315bcSPatrick Sanan             /* left-hand top down corner */
690e77315bcSPatrick Sanan             if (k == 0) {
691e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
692e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
693e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
694e77315bcSPatrick Sanan               /* du = bu * au; */
695e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
696e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
697e77315bcSPatrick Sanan 
6989371c9d4SSatish Balay               c[0].k = k;
6999371c9d4SSatish Balay               c[0].j = j - 1;
7009371c9d4SSatish Balay               c[0].i = i;
701e77315bcSPatrick Sanan               v[0]   = -hzhxdhy * (ds - gs);
7029371c9d4SSatish Balay               c[1].k = k;
7039371c9d4SSatish Balay               c[1].j = j;
7049371c9d4SSatish Balay               c[1].i = i;
705e77315bcSPatrick Sanan               v[1]   = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu);
7069371c9d4SSatish Balay               c[2].k = k;
7079371c9d4SSatish Balay               c[2].j = j;
7089371c9d4SSatish Balay               c[2].i = i + 1;
709e77315bcSPatrick Sanan               v[2]   = -hyhzdhx * (de + ge);
7109371c9d4SSatish Balay               c[3].k = k + 1;
7119371c9d4SSatish Balay               c[3].j = j;
7129371c9d4SSatish Balay               c[3].i = i;
713e77315bcSPatrick Sanan               v[3]   = -hxhydhz * (du + gu);
7149566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES));
715e77315bcSPatrick Sanan 
716e77315bcSPatrick Sanan               /* left-hand top interior edge */
717e77315bcSPatrick Sanan             } else if (k < mz - 1) {
718e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
719e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
720e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
721e77315bcSPatrick Sanan               /* du = bu * au; */
722e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
723e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
724e77315bcSPatrick Sanan 
725e77315bcSPatrick Sanan               td = x[k - 1][j][i];
726e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
727e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
728e77315bcSPatrick Sanan               /* dd = bd * ad; */
729e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
730e77315bcSPatrick Sanan               gd = coef * bd * (td - t0);
731e77315bcSPatrick Sanan 
7329371c9d4SSatish Balay               c[0].k = k - 1;
7339371c9d4SSatish Balay               c[0].j = j;
7349371c9d4SSatish Balay               c[0].i = i;
735e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
7369371c9d4SSatish Balay               c[1].k = k;
7379371c9d4SSatish Balay               c[1].j = j - 1;
7389371c9d4SSatish Balay               c[1].i = i;
739e77315bcSPatrick Sanan               v[1]   = -hzhxdhy * (ds - gs);
7409371c9d4SSatish Balay               c[2].k = k;
7419371c9d4SSatish Balay               c[2].j = j;
7429371c9d4SSatish Balay               c[2].i = i;
743e77315bcSPatrick Sanan               v[2]   = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu);
7449371c9d4SSatish Balay               c[3].k = k;
7459371c9d4SSatish Balay               c[3].j = j;
7469371c9d4SSatish Balay               c[3].i = i + 1;
747e77315bcSPatrick Sanan               v[3]   = -hyhzdhx * (de + ge);
7489371c9d4SSatish Balay               c[4].k = k + 1;
7499371c9d4SSatish Balay               c[4].j = j;
7509371c9d4SSatish Balay               c[4].i = i;
751e77315bcSPatrick Sanan               v[4]   = -hxhydhz * (du + gu);
7529566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
753e77315bcSPatrick Sanan 
754e77315bcSPatrick Sanan               /* left-hand top up corner */
755e77315bcSPatrick Sanan             } else {
756e77315bcSPatrick Sanan               td = x[k - 1][j][i];
757e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
758e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
759e77315bcSPatrick Sanan               /* dd = bd * ad; */
760e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
761e77315bcSPatrick Sanan               gd = coef * bd * (td - t0);
762e77315bcSPatrick Sanan 
7639371c9d4SSatish Balay               c[0].k = k - 1;
7649371c9d4SSatish Balay               c[0].j = j;
7659371c9d4SSatish Balay               c[0].i = i;
766e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
7679371c9d4SSatish Balay               c[1].k = k;
7689371c9d4SSatish Balay               c[1].j = j - 1;
7699371c9d4SSatish Balay               c[1].i = i;
770e77315bcSPatrick Sanan               v[1]   = -hzhxdhy * (ds - gs);
7719371c9d4SSatish Balay               c[2].k = k;
7729371c9d4SSatish Balay               c[2].j = j;
7739371c9d4SSatish Balay               c[2].i = i;
774e77315bcSPatrick Sanan               v[2]   = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd);
7759371c9d4SSatish Balay               c[3].k = k;
7769371c9d4SSatish Balay               c[3].j = j;
7779371c9d4SSatish Balay               c[3].i = i + 1;
778e77315bcSPatrick Sanan               v[3]   = -hyhzdhx * (de + ge);
7799566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES));
780e77315bcSPatrick Sanan             }
781e77315bcSPatrick Sanan 
782e77315bcSPatrick Sanan           } else {
783e77315bcSPatrick Sanan             ts = x[k][j - 1][i];
784e77315bcSPatrick Sanan             as = 0.5 * (t0 + ts);
785e77315bcSPatrick Sanan             bs = PetscPowScalar(as, bm1);
786e77315bcSPatrick Sanan             /* ds = bs * as; */
787e77315bcSPatrick Sanan             ds = PetscPowScalar(as, beta);
788e77315bcSPatrick Sanan             gs = coef * bs * (t0 - ts);
789e77315bcSPatrick Sanan 
790e77315bcSPatrick Sanan             tn = x[k][j + 1][i];
791e77315bcSPatrick Sanan             an = 0.5 * (t0 + tn);
792e77315bcSPatrick Sanan             bn = PetscPowScalar(an, bm1);
793e77315bcSPatrick Sanan             /* dn = bn * an; */
794e77315bcSPatrick Sanan             dn = PetscPowScalar(an, beta);
795e77315bcSPatrick Sanan             gn = coef * bn * (tn - t0);
796e77315bcSPatrick Sanan 
797e77315bcSPatrick Sanan             /* left-hand down interior edge */
798e77315bcSPatrick Sanan             if (k == 0) {
799e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
800e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
801e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
802e77315bcSPatrick Sanan               /* du = bu * au; */
803e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
804e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
805e77315bcSPatrick Sanan 
8069371c9d4SSatish Balay               c[0].k = k;
8079371c9d4SSatish Balay               c[0].j = j - 1;
8089371c9d4SSatish Balay               c[0].i = i;
809e77315bcSPatrick Sanan               v[0]   = -hzhxdhy * (ds - gs);
8109371c9d4SSatish Balay               c[1].k = k;
8119371c9d4SSatish Balay               c[1].j = j;
8129371c9d4SSatish Balay               c[1].i = i;
813e77315bcSPatrick Sanan               v[1]   = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu);
8149371c9d4SSatish Balay               c[2].k = k;
8159371c9d4SSatish Balay               c[2].j = j;
8169371c9d4SSatish Balay               c[2].i = i + 1;
817e77315bcSPatrick Sanan               v[2]   = -hyhzdhx * (de + ge);
8189371c9d4SSatish Balay               c[3].k = k;
8199371c9d4SSatish Balay               c[3].j = j + 1;
8209371c9d4SSatish Balay               c[3].i = i;
821e77315bcSPatrick Sanan               v[3]   = -hzhxdhy * (dn + gn);
8229371c9d4SSatish Balay               c[4].k = k + 1;
8239371c9d4SSatish Balay               c[4].j = j;
8249371c9d4SSatish Balay               c[4].i = i;
825e77315bcSPatrick Sanan               v[4]   = -hxhydhz * (du + gu);
8269566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
827e77315bcSPatrick Sanan 
828e77315bcSPatrick Sanan             } else if (k == mz - 1) { /* left-hand up interior edge */
829e77315bcSPatrick Sanan 
830e77315bcSPatrick Sanan               td = x[k - 1][j][i];
831e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
832e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
833e77315bcSPatrick Sanan               /* dd = bd * ad; */
834e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
835e77315bcSPatrick Sanan               gd = coef * bd * (t0 - td);
836e77315bcSPatrick Sanan 
8379371c9d4SSatish Balay               c[0].k = k - 1;
8389371c9d4SSatish Balay               c[0].j = j;
8399371c9d4SSatish Balay               c[0].i = i;
840e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
8419371c9d4SSatish Balay               c[1].k = k;
8429371c9d4SSatish Balay               c[1].j = j - 1;
8439371c9d4SSatish Balay               c[1].i = i;
844e77315bcSPatrick Sanan               v[1]   = -hzhxdhy * (ds - gs);
8459371c9d4SSatish Balay               c[2].k = k;
8469371c9d4SSatish Balay               c[2].j = j;
8479371c9d4SSatish Balay               c[2].i = i;
848e77315bcSPatrick Sanan               v[2]   = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd);
8499371c9d4SSatish Balay               c[3].k = k;
8509371c9d4SSatish Balay               c[3].j = j;
8519371c9d4SSatish Balay               c[3].i = i + 1;
852e77315bcSPatrick Sanan               v[3]   = -hyhzdhx * (de + ge);
8539371c9d4SSatish Balay               c[4].k = k;
8549371c9d4SSatish Balay               c[4].j = j + 1;
8559371c9d4SSatish Balay               c[4].i = i;
856e77315bcSPatrick Sanan               v[4]   = -hzhxdhy * (dn + gn);
8579566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
858e77315bcSPatrick Sanan             } else { /* left-hand interior plane */
859e77315bcSPatrick Sanan 
860e77315bcSPatrick Sanan               td = x[k - 1][j][i];
861e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
862e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
863e77315bcSPatrick Sanan               /* dd = bd * ad; */
864e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
865e77315bcSPatrick Sanan               gd = coef * bd * (t0 - td);
866e77315bcSPatrick Sanan 
867e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
868e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
869e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
870e77315bcSPatrick Sanan               /* du = bu * au; */
871e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
872e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
873e77315bcSPatrick Sanan 
8749371c9d4SSatish Balay               c[0].k = k - 1;
8759371c9d4SSatish Balay               c[0].j = j;
8769371c9d4SSatish Balay               c[0].i = i;
877e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
8789371c9d4SSatish Balay               c[1].k = k;
8799371c9d4SSatish Balay               c[1].j = j - 1;
8809371c9d4SSatish Balay               c[1].i = i;
881e77315bcSPatrick Sanan               v[1]   = -hzhxdhy * (ds - gs);
8829371c9d4SSatish Balay               c[2].k = k;
8839371c9d4SSatish Balay               c[2].j = j;
8849371c9d4SSatish Balay               c[2].i = i;
885e77315bcSPatrick Sanan               v[2]   = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu);
8869371c9d4SSatish Balay               c[3].k = k;
8879371c9d4SSatish Balay               c[3].j = j;
8889371c9d4SSatish Balay               c[3].i = i + 1;
889e77315bcSPatrick Sanan               v[3]   = -hyhzdhx * (de + ge);
8909371c9d4SSatish Balay               c[4].k = k;
8919371c9d4SSatish Balay               c[4].j = j + 1;
8929371c9d4SSatish Balay               c[4].i = i;
893e77315bcSPatrick Sanan               v[4]   = -hzhxdhy * (dn + gn);
8949371c9d4SSatish Balay               c[5].k = k + 1;
8959371c9d4SSatish Balay               c[5].j = j;
8969371c9d4SSatish Balay               c[5].i = i;
897e77315bcSPatrick Sanan               v[5]   = -hxhydhz * (du + gu);
8989566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES));
899e77315bcSPatrick Sanan             }
900e77315bcSPatrick Sanan           }
901e77315bcSPatrick Sanan 
902e77315bcSPatrick Sanan         } else if (i == mx - 1) {
903e77315bcSPatrick Sanan           /* right-hand plane boundary */
904e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
905e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
906e77315bcSPatrick Sanan           bw = PetscPowScalar(aw, bm1);
907e77315bcSPatrick Sanan           /* dw = bw * aw */
908e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
909e77315bcSPatrick Sanan           gw = coef * bw * (t0 - tw);
910e77315bcSPatrick Sanan 
911e77315bcSPatrick Sanan           te = tright;
912e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
913e77315bcSPatrick Sanan           be = PetscPowScalar(ae, bm1);
914e77315bcSPatrick Sanan           /* de = be * ae; */
915e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
916e77315bcSPatrick Sanan           ge = coef * be * (te - t0);
917e77315bcSPatrick Sanan 
918e77315bcSPatrick Sanan           /* right-hand bottom edge */
919e77315bcSPatrick Sanan           if (j == 0) {
920e77315bcSPatrick Sanan             tn = x[k][j + 1][i];
921e77315bcSPatrick Sanan             an = 0.5 * (t0 + tn);
922e77315bcSPatrick Sanan             bn = PetscPowScalar(an, bm1);
923e77315bcSPatrick Sanan             /* dn = bn * an; */
924e77315bcSPatrick Sanan             dn = PetscPowScalar(an, beta);
925e77315bcSPatrick Sanan             gn = coef * bn * (tn - t0);
926e77315bcSPatrick Sanan 
927e77315bcSPatrick Sanan             /* right-hand bottom down corner */
928e77315bcSPatrick Sanan             if (k == 0) {
929e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
930e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
931e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
932e77315bcSPatrick Sanan               /* du = bu * au; */
933e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
934e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
935e77315bcSPatrick Sanan 
9369371c9d4SSatish Balay               c[0].k = k;
9379371c9d4SSatish Balay               c[0].j = j;
9389371c9d4SSatish Balay               c[0].i = i - 1;
939e77315bcSPatrick Sanan               v[0]   = -hyhzdhx * (dw - gw);
9409371c9d4SSatish Balay               c[1].k = k;
9419371c9d4SSatish Balay               c[1].j = j;
9429371c9d4SSatish Balay               c[1].i = i;
943e77315bcSPatrick Sanan               v[1]   = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu);
9449371c9d4SSatish Balay               c[2].k = k;
9459371c9d4SSatish Balay               c[2].j = j + 1;
9469371c9d4SSatish Balay               c[2].i = i;
947e77315bcSPatrick Sanan               v[2]   = -hzhxdhy * (dn + gn);
9489371c9d4SSatish Balay               c[3].k = k + 1;
9499371c9d4SSatish Balay               c[3].j = j;
9509371c9d4SSatish Balay               c[3].i = i;
951e77315bcSPatrick Sanan               v[3]   = -hxhydhz * (du + gu);
9529566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES));
953e77315bcSPatrick Sanan 
954e77315bcSPatrick Sanan               /* right-hand bottom interior edge */
955e77315bcSPatrick Sanan             } else if (k < mz - 1) {
956e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
957e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
958e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
959e77315bcSPatrick Sanan               /* du = bu * au; */
960e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
961e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
962e77315bcSPatrick Sanan 
963e77315bcSPatrick Sanan               td = x[k - 1][j][i];
964e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
965e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
966e77315bcSPatrick Sanan               /* dd = bd * ad; */
967e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
968e77315bcSPatrick Sanan               gd = coef * bd * (td - t0);
969e77315bcSPatrick Sanan 
9709371c9d4SSatish Balay               c[0].k = k - 1;
9719371c9d4SSatish Balay               c[0].j = j;
9729371c9d4SSatish Balay               c[0].i = i;
973e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
9749371c9d4SSatish Balay               c[1].k = k;
9759371c9d4SSatish Balay               c[1].j = j;
9769371c9d4SSatish Balay               c[1].i = i - 1;
977e77315bcSPatrick Sanan               v[1]   = -hyhzdhx * (dw - gw);
9789371c9d4SSatish Balay               c[2].k = k;
9799371c9d4SSatish Balay               c[2].j = j;
9809371c9d4SSatish Balay               c[2].i = i;
981e77315bcSPatrick Sanan               v[2]   = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu);
9829371c9d4SSatish Balay               c[3].k = k;
9839371c9d4SSatish Balay               c[3].j = j + 1;
9849371c9d4SSatish Balay               c[3].i = i;
985e77315bcSPatrick Sanan               v[3]   = -hzhxdhy * (dn + gn);
9869371c9d4SSatish Balay               c[4].k = k + 1;
9879371c9d4SSatish Balay               c[4].j = j;
9889371c9d4SSatish Balay               c[4].i = i;
989e77315bcSPatrick Sanan               v[4]   = -hxhydhz * (du + gu);
9909566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
991e77315bcSPatrick Sanan 
992e77315bcSPatrick Sanan               /* right-hand bottom up corner */
993e77315bcSPatrick Sanan             } else {
994e77315bcSPatrick Sanan               td = x[k - 1][j][i];
995e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
996e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
997e77315bcSPatrick Sanan               /* dd = bd * ad; */
998e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
999e77315bcSPatrick Sanan               gd = coef * bd * (td - t0);
1000e77315bcSPatrick Sanan 
10019371c9d4SSatish Balay               c[0].k = k - 1;
10029371c9d4SSatish Balay               c[0].j = j;
10039371c9d4SSatish Balay               c[0].i = i;
1004e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
10059371c9d4SSatish Balay               c[1].k = k;
10069371c9d4SSatish Balay               c[1].j = j;
10079371c9d4SSatish Balay               c[1].i = i - 1;
1008e77315bcSPatrick Sanan               v[1]   = -hyhzdhx * (dw - gw);
10099371c9d4SSatish Balay               c[2].k = k;
10109371c9d4SSatish Balay               c[2].j = j;
10119371c9d4SSatish Balay               c[2].i = i;
1012e77315bcSPatrick Sanan               v[2]   = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd);
10139371c9d4SSatish Balay               c[3].k = k;
10149371c9d4SSatish Balay               c[3].j = j + 1;
10159371c9d4SSatish Balay               c[3].i = i;
1016e77315bcSPatrick Sanan               v[3]   = -hzhxdhy * (dn + gn);
10179566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES));
1018e77315bcSPatrick Sanan             }
1019e77315bcSPatrick Sanan 
1020e77315bcSPatrick Sanan             /* right-hand top edge */
1021e77315bcSPatrick Sanan           } else if (j == my - 1) {
1022e77315bcSPatrick Sanan             ts = x[k][j - 1][i];
1023e77315bcSPatrick Sanan             as = 0.5 * (t0 + ts);
1024e77315bcSPatrick Sanan             bs = PetscPowScalar(as, bm1);
1025e77315bcSPatrick Sanan             /* ds = bs * as; */
1026e77315bcSPatrick Sanan             ds = PetscPowScalar(as, beta);
1027e77315bcSPatrick Sanan             gs = coef * bs * (ts - t0);
1028e77315bcSPatrick Sanan 
1029e77315bcSPatrick Sanan             /* right-hand top down corner */
1030e77315bcSPatrick Sanan             if (k == 0) {
1031e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
1032e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
1033e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
1034e77315bcSPatrick Sanan               /* du = bu * au; */
1035e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
1036e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
1037e77315bcSPatrick Sanan 
10389371c9d4SSatish Balay               c[0].k = k;
10399371c9d4SSatish Balay               c[0].j = j - 1;
10409371c9d4SSatish Balay               c[0].i = i;
1041e77315bcSPatrick Sanan               v[0]   = -hzhxdhy * (ds - gs);
10429371c9d4SSatish Balay               c[1].k = k;
10439371c9d4SSatish Balay               c[1].j = j;
10449371c9d4SSatish Balay               c[1].i = i - 1;
1045e77315bcSPatrick Sanan               v[1]   = -hyhzdhx * (dw - gw);
10469371c9d4SSatish Balay               c[2].k = k;
10479371c9d4SSatish Balay               c[2].j = j;
10489371c9d4SSatish Balay               c[2].i = i;
1049e77315bcSPatrick Sanan               v[2]   = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu);
10509371c9d4SSatish Balay               c[3].k = k + 1;
10519371c9d4SSatish Balay               c[3].j = j;
10529371c9d4SSatish Balay               c[3].i = i;
1053e77315bcSPatrick Sanan               v[3]   = -hxhydhz * (du + gu);
10549566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES));
1055e77315bcSPatrick Sanan 
1056e77315bcSPatrick Sanan               /* right-hand top interior edge */
1057e77315bcSPatrick Sanan             } else if (k < mz - 1) {
1058e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
1059e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
1060e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
1061e77315bcSPatrick Sanan               /* du = bu * au; */
1062e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
1063e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
1064e77315bcSPatrick Sanan 
1065e77315bcSPatrick Sanan               td = x[k - 1][j][i];
1066e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
1067e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
1068e77315bcSPatrick Sanan               /* dd = bd * ad; */
1069e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
1070e77315bcSPatrick Sanan               gd = coef * bd * (td - t0);
1071e77315bcSPatrick Sanan 
10729371c9d4SSatish Balay               c[0].k = k - 1;
10739371c9d4SSatish Balay               c[0].j = j;
10749371c9d4SSatish Balay               c[0].i = i;
1075e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
10769371c9d4SSatish Balay               c[1].k = k;
10779371c9d4SSatish Balay               c[1].j = j - 1;
10789371c9d4SSatish Balay               c[1].i = i;
1079e77315bcSPatrick Sanan               v[1]   = -hzhxdhy * (ds - gs);
10809371c9d4SSatish Balay               c[2].k = k;
10819371c9d4SSatish Balay               c[2].j = j;
10829371c9d4SSatish Balay               c[2].i = i - 1;
1083e77315bcSPatrick Sanan               v[2]   = -hyhzdhx * (dw - gw);
10849371c9d4SSatish Balay               c[3].k = k;
10859371c9d4SSatish Balay               c[3].j = j;
10869371c9d4SSatish Balay               c[3].i = i;
1087e77315bcSPatrick Sanan               v[3]   = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu);
10889371c9d4SSatish Balay               c[4].k = k + 1;
10899371c9d4SSatish Balay               c[4].j = j;
10909371c9d4SSatish Balay               c[4].i = i;
1091e77315bcSPatrick Sanan               v[4]   = -hxhydhz * (du + gu);
10929566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
1093e77315bcSPatrick Sanan 
1094e77315bcSPatrick Sanan               /* right-hand top up corner */
1095e77315bcSPatrick Sanan             } else {
1096e77315bcSPatrick Sanan               td = x[k - 1][j][i];
1097e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
1098e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
1099e77315bcSPatrick Sanan               /* dd = bd * ad; */
1100e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
1101e77315bcSPatrick Sanan               gd = coef * bd * (td - t0);
1102e77315bcSPatrick Sanan 
11039371c9d4SSatish Balay               c[0].k = k - 1;
11049371c9d4SSatish Balay               c[0].j = j;
11059371c9d4SSatish Balay               c[0].i = i;
1106e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
11079371c9d4SSatish Balay               c[1].k = k;
11089371c9d4SSatish Balay               c[1].j = j - 1;
11099371c9d4SSatish Balay               c[1].i = i;
1110e77315bcSPatrick Sanan               v[1]   = -hzhxdhy * (ds - gs);
11119371c9d4SSatish Balay               c[2].k = k;
11129371c9d4SSatish Balay               c[2].j = j;
11139371c9d4SSatish Balay               c[2].i = i - 1;
1114e77315bcSPatrick Sanan               v[2]   = -hyhzdhx * (dw - gw);
11159371c9d4SSatish Balay               c[3].k = k;
11169371c9d4SSatish Balay               c[3].j = j;
11179371c9d4SSatish Balay               c[3].i = i;
1118e77315bcSPatrick Sanan               v[3]   = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd);
11199566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 4, c, v, INSERT_VALUES));
1120e77315bcSPatrick Sanan             }
1121e77315bcSPatrick Sanan 
1122e77315bcSPatrick Sanan           } else {
1123e77315bcSPatrick Sanan             ts = x[k][j - 1][i];
1124e77315bcSPatrick Sanan             as = 0.5 * (t0 + ts);
1125e77315bcSPatrick Sanan             bs = PetscPowScalar(as, bm1);
1126e77315bcSPatrick Sanan             /* ds = bs * as; */
1127e77315bcSPatrick Sanan             ds = PetscPowScalar(as, beta);
1128e77315bcSPatrick Sanan             gs = coef * bs * (t0 - ts);
1129e77315bcSPatrick Sanan 
1130e77315bcSPatrick Sanan             tn = x[k][j + 1][i];
1131e77315bcSPatrick Sanan             an = 0.5 * (t0 + tn);
1132e77315bcSPatrick Sanan             bn = PetscPowScalar(an, bm1);
1133e77315bcSPatrick Sanan             /* dn = bn * an; */
1134e77315bcSPatrick Sanan             dn = PetscPowScalar(an, beta);
1135e77315bcSPatrick Sanan             gn = coef * bn * (tn - t0);
1136e77315bcSPatrick Sanan 
1137e77315bcSPatrick Sanan             /* right-hand down interior edge */
1138e77315bcSPatrick Sanan             if (k == 0) {
1139e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
1140e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
1141e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
1142e77315bcSPatrick Sanan               /* du = bu * au; */
1143e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
1144e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
1145e77315bcSPatrick Sanan 
11469371c9d4SSatish Balay               c[0].k = k;
11479371c9d4SSatish Balay               c[0].j = j - 1;
11489371c9d4SSatish Balay               c[0].i = i;
1149e77315bcSPatrick Sanan               v[0]   = -hzhxdhy * (ds - gs);
11509371c9d4SSatish Balay               c[1].k = k;
11519371c9d4SSatish Balay               c[1].j = j;
11529371c9d4SSatish Balay               c[1].i = i - 1;
1153e77315bcSPatrick Sanan               v[1]   = -hyhzdhx * (dw - gw);
11549371c9d4SSatish Balay               c[2].k = k;
11559371c9d4SSatish Balay               c[2].j = j;
11569371c9d4SSatish Balay               c[2].i = i;
1157e77315bcSPatrick Sanan               v[2]   = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu);
11589371c9d4SSatish Balay               c[3].k = k;
11599371c9d4SSatish Balay               c[3].j = j + 1;
11609371c9d4SSatish Balay               c[3].i = i;
1161e77315bcSPatrick Sanan               v[3]   = -hzhxdhy * (dn + gn);
11629371c9d4SSatish Balay               c[4].k = k + 1;
11639371c9d4SSatish Balay               c[4].j = j;
11649371c9d4SSatish Balay               c[4].i = i;
1165e77315bcSPatrick Sanan               v[4]   = -hxhydhz * (du + gu);
11669566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
1167e77315bcSPatrick Sanan 
1168e77315bcSPatrick Sanan             } else if (k == mz - 1) { /* right-hand up interior edge */
1169e77315bcSPatrick Sanan 
1170e77315bcSPatrick Sanan               td = x[k - 1][j][i];
1171e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
1172e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
1173e77315bcSPatrick Sanan               /* dd = bd * ad; */
1174e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
1175e77315bcSPatrick Sanan               gd = coef * bd * (t0 - td);
1176e77315bcSPatrick Sanan 
11779371c9d4SSatish Balay               c[0].k = k - 1;
11789371c9d4SSatish Balay               c[0].j = j;
11799371c9d4SSatish Balay               c[0].i = i;
1180e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
11819371c9d4SSatish Balay               c[1].k = k;
11829371c9d4SSatish Balay               c[1].j = j - 1;
11839371c9d4SSatish Balay               c[1].i = i;
1184e77315bcSPatrick Sanan               v[1]   = -hzhxdhy * (ds - gs);
11859371c9d4SSatish Balay               c[2].k = k;
11869371c9d4SSatish Balay               c[2].j = j;
11879371c9d4SSatish Balay               c[2].i = i - 1;
1188e77315bcSPatrick Sanan               v[2]   = -hyhzdhx * (dw - gw);
11899371c9d4SSatish Balay               c[3].k = k;
11909371c9d4SSatish Balay               c[3].j = j;
11919371c9d4SSatish Balay               c[3].i = i;
1192e77315bcSPatrick Sanan               v[3]   = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd);
11939371c9d4SSatish Balay               c[4].k = k;
11949371c9d4SSatish Balay               c[4].j = j + 1;
11959371c9d4SSatish Balay               c[4].i = i;
1196e77315bcSPatrick Sanan               v[4]   = -hzhxdhy * (dn + gn);
11979566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
1198e77315bcSPatrick Sanan 
1199e77315bcSPatrick Sanan             } else { /* right-hand interior plane */
1200e77315bcSPatrick Sanan 
1201e77315bcSPatrick Sanan               td = x[k - 1][j][i];
1202e77315bcSPatrick Sanan               ad = 0.5 * (t0 + td);
1203e77315bcSPatrick Sanan               bd = PetscPowScalar(ad, bm1);
1204e77315bcSPatrick Sanan               /* dd = bd * ad; */
1205e77315bcSPatrick Sanan               dd = PetscPowScalar(ad, beta);
1206e77315bcSPatrick Sanan               gd = coef * bd * (t0 - td);
1207e77315bcSPatrick Sanan 
1208e77315bcSPatrick Sanan               tu = x[k + 1][j][i];
1209e77315bcSPatrick Sanan               au = 0.5 * (t0 + tu);
1210e77315bcSPatrick Sanan               bu = PetscPowScalar(au, bm1);
1211e77315bcSPatrick Sanan               /* du = bu * au; */
1212e77315bcSPatrick Sanan               du = PetscPowScalar(au, beta);
1213e77315bcSPatrick Sanan               gu = coef * bu * (tu - t0);
1214e77315bcSPatrick Sanan 
12159371c9d4SSatish Balay               c[0].k = k - 1;
12169371c9d4SSatish Balay               c[0].j = j;
12179371c9d4SSatish Balay               c[0].i = i;
1218e77315bcSPatrick Sanan               v[0]   = -hxhydhz * (dd - gd);
12199371c9d4SSatish Balay               c[1].k = k;
12209371c9d4SSatish Balay               c[1].j = j - 1;
12219371c9d4SSatish Balay               c[1].i = i;
1222e77315bcSPatrick Sanan               v[1]   = -hzhxdhy * (ds - gs);
12239371c9d4SSatish Balay               c[2].k = k;
12249371c9d4SSatish Balay               c[2].j = j;
12259371c9d4SSatish Balay               c[2].i = i - 1;
1226e77315bcSPatrick Sanan               v[2]   = -hyhzdhx * (dw - gw);
12279371c9d4SSatish Balay               c[3].k = k;
12289371c9d4SSatish Balay               c[3].j = j;
12299371c9d4SSatish Balay               c[3].i = i;
1230e77315bcSPatrick Sanan               v[3]   = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu);
12319371c9d4SSatish Balay               c[4].k = k;
12329371c9d4SSatish Balay               c[4].j = j + 1;
12339371c9d4SSatish Balay               c[4].i = i;
1234e77315bcSPatrick Sanan               v[4]   = -hzhxdhy * (dn + gn);
12359371c9d4SSatish Balay               c[5].k = k + 1;
12369371c9d4SSatish Balay               c[5].j = j;
12379371c9d4SSatish Balay               c[5].i = i;
1238e77315bcSPatrick Sanan               v[5]   = -hxhydhz * (du + gu);
12399566063dSJacob Faibussowitsch               PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES));
1240e77315bcSPatrick Sanan             }
1241e77315bcSPatrick Sanan           }
1242e77315bcSPatrick Sanan 
1243e77315bcSPatrick Sanan         } else if (j == 0) {
1244e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
1245e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
1246e77315bcSPatrick Sanan           bw = PetscPowScalar(aw, bm1);
1247e77315bcSPatrick Sanan           /* dw = bw * aw */
1248e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
1249e77315bcSPatrick Sanan           gw = coef * bw * (t0 - tw);
1250e77315bcSPatrick Sanan 
1251e77315bcSPatrick Sanan           te = x[k][j][i + 1];
1252e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
1253e77315bcSPatrick Sanan           be = PetscPowScalar(ae, bm1);
1254e77315bcSPatrick Sanan           /* de = be * ae; */
1255e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
1256e77315bcSPatrick Sanan           ge = coef * be * (te - t0);
1257e77315bcSPatrick Sanan 
1258e77315bcSPatrick Sanan           tn = x[k][j + 1][i];
1259e77315bcSPatrick Sanan           an = 0.5 * (t0 + tn);
1260e77315bcSPatrick Sanan           bn = PetscPowScalar(an, bm1);
1261e77315bcSPatrick Sanan           /* dn = bn * an; */
1262e77315bcSPatrick Sanan           dn = PetscPowScalar(an, beta);
1263e77315bcSPatrick Sanan           gn = coef * bn * (tn - t0);
1264e77315bcSPatrick Sanan 
1265e77315bcSPatrick Sanan           /* bottom down interior edge */
1266e77315bcSPatrick Sanan           if (k == 0) {
1267e77315bcSPatrick Sanan             tu = x[k + 1][j][i];
1268e77315bcSPatrick Sanan             au = 0.5 * (t0 + tu);
1269e77315bcSPatrick Sanan             bu = PetscPowScalar(au, bm1);
1270e77315bcSPatrick Sanan             /* du = bu * au; */
1271e77315bcSPatrick Sanan             du = PetscPowScalar(au, beta);
1272e77315bcSPatrick Sanan             gu = coef * bu * (tu - t0);
1273e77315bcSPatrick Sanan 
12749371c9d4SSatish Balay             c[0].k = k;
12759371c9d4SSatish Balay             c[0].j = j;
12769371c9d4SSatish Balay             c[0].i = i - 1;
1277e77315bcSPatrick Sanan             v[0]   = -hyhzdhx * (dw - gw);
12789371c9d4SSatish Balay             c[1].k = k;
12799371c9d4SSatish Balay             c[1].j = j;
12809371c9d4SSatish Balay             c[1].i = i;
1281e77315bcSPatrick Sanan             v[1]   = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu);
12829371c9d4SSatish Balay             c[2].k = k;
12839371c9d4SSatish Balay             c[2].j = j;
12849371c9d4SSatish Balay             c[2].i = i + 1;
1285e77315bcSPatrick Sanan             v[2]   = -hyhzdhx * (de + ge);
12869371c9d4SSatish Balay             c[3].k = k;
12879371c9d4SSatish Balay             c[3].j = j + 1;
12889371c9d4SSatish Balay             c[3].i = i;
1289e77315bcSPatrick Sanan             v[3]   = -hzhxdhy * (dn + gn);
12909371c9d4SSatish Balay             c[4].k = k + 1;
12919371c9d4SSatish Balay             c[4].j = j;
12929371c9d4SSatish Balay             c[4].i = i;
1293e77315bcSPatrick Sanan             v[4]   = -hxhydhz * (du + gu);
12949566063dSJacob Faibussowitsch             PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
1295e77315bcSPatrick Sanan 
1296e77315bcSPatrick Sanan           } else if (k == mz - 1) { /* bottom up interior edge */
1297e77315bcSPatrick Sanan 
1298e77315bcSPatrick Sanan             td = x[k - 1][j][i];
1299e77315bcSPatrick Sanan             ad = 0.5 * (t0 + td);
1300e77315bcSPatrick Sanan             bd = PetscPowScalar(ad, bm1);
1301e77315bcSPatrick Sanan             /* dd = bd * ad; */
1302e77315bcSPatrick Sanan             dd = PetscPowScalar(ad, beta);
1303e77315bcSPatrick Sanan             gd = coef * bd * (td - t0);
1304e77315bcSPatrick Sanan 
13059371c9d4SSatish Balay             c[0].k = k - 1;
13069371c9d4SSatish Balay             c[0].j = j;
13079371c9d4SSatish Balay             c[0].i = i;
1308e77315bcSPatrick Sanan             v[0]   = -hxhydhz * (dd - gd);
13099371c9d4SSatish Balay             c[1].k = k;
13109371c9d4SSatish Balay             c[1].j = j;
13119371c9d4SSatish Balay             c[1].i = i - 1;
1312e77315bcSPatrick Sanan             v[1]   = -hyhzdhx * (dw - gw);
13139371c9d4SSatish Balay             c[2].k = k;
13149371c9d4SSatish Balay             c[2].j = j;
13159371c9d4SSatish Balay             c[2].i = i;
1316e77315bcSPatrick Sanan             v[2]   = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd);
13179371c9d4SSatish Balay             c[3].k = k;
13189371c9d4SSatish Balay             c[3].j = j;
13199371c9d4SSatish Balay             c[3].i = i + 1;
1320e77315bcSPatrick Sanan             v[3]   = -hyhzdhx * (de + ge);
13219371c9d4SSatish Balay             c[4].k = k;
13229371c9d4SSatish Balay             c[4].j = j + 1;
13239371c9d4SSatish Balay             c[4].i = i;
1324e77315bcSPatrick Sanan             v[4]   = -hzhxdhy * (dn + gn);
13259566063dSJacob Faibussowitsch             PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
1326e77315bcSPatrick Sanan 
1327e77315bcSPatrick Sanan           } else { /* bottom interior plane */
1328e77315bcSPatrick Sanan 
1329e77315bcSPatrick Sanan             tu = x[k + 1][j][i];
1330e77315bcSPatrick Sanan             au = 0.5 * (t0 + tu);
1331e77315bcSPatrick Sanan             bu = PetscPowScalar(au, bm1);
1332e77315bcSPatrick Sanan             /* du = bu * au; */
1333e77315bcSPatrick Sanan             du = PetscPowScalar(au, beta);
1334e77315bcSPatrick Sanan             gu = coef * bu * (tu - t0);
1335e77315bcSPatrick Sanan 
1336e77315bcSPatrick Sanan             td = x[k - 1][j][i];
1337e77315bcSPatrick Sanan             ad = 0.5 * (t0 + td);
1338e77315bcSPatrick Sanan             bd = PetscPowScalar(ad, bm1);
1339e77315bcSPatrick Sanan             /* dd = bd * ad; */
1340e77315bcSPatrick Sanan             dd = PetscPowScalar(ad, beta);
1341e77315bcSPatrick Sanan             gd = coef * bd * (td - t0);
1342e77315bcSPatrick Sanan 
13439371c9d4SSatish Balay             c[0].k = k - 1;
13449371c9d4SSatish Balay             c[0].j = j;
13459371c9d4SSatish Balay             c[0].i = i;
1346e77315bcSPatrick Sanan             v[0]   = -hxhydhz * (dd - gd);
13479371c9d4SSatish Balay             c[1].k = k;
13489371c9d4SSatish Balay             c[1].j = j;
13499371c9d4SSatish Balay             c[1].i = i - 1;
1350e77315bcSPatrick Sanan             v[1]   = -hyhzdhx * (dw - gw);
13519371c9d4SSatish Balay             c[2].k = k;
13529371c9d4SSatish Balay             c[2].j = j;
13539371c9d4SSatish Balay             c[2].i = i;
1354e77315bcSPatrick Sanan             v[2]   = hzhxdhy * (dn - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu);
13559371c9d4SSatish Balay             c[3].k = k;
13569371c9d4SSatish Balay             c[3].j = j;
13579371c9d4SSatish Balay             c[3].i = i + 1;
1358e77315bcSPatrick Sanan             v[3]   = -hyhzdhx * (de + ge);
13599371c9d4SSatish Balay             c[4].k = k;
13609371c9d4SSatish Balay             c[4].j = j + 1;
13619371c9d4SSatish Balay             c[4].i = i;
1362e77315bcSPatrick Sanan             v[4]   = -hzhxdhy * (dn + gn);
13639371c9d4SSatish Balay             c[5].k = k + 1;
13649371c9d4SSatish Balay             c[5].j = j;
13659371c9d4SSatish Balay             c[5].i = i;
1366e77315bcSPatrick Sanan             v[5]   = -hxhydhz * (du + gu);
13679566063dSJacob Faibussowitsch             PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES));
1368e77315bcSPatrick Sanan           }
1369e77315bcSPatrick Sanan 
1370e77315bcSPatrick Sanan         } else if (j == my - 1) {
1371e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
1372e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
1373e77315bcSPatrick Sanan           bw = PetscPowScalar(aw, bm1);
1374e77315bcSPatrick Sanan           /* dw = bw * aw */
1375e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
1376e77315bcSPatrick Sanan           gw = coef * bw * (t0 - tw);
1377e77315bcSPatrick Sanan 
1378e77315bcSPatrick Sanan           te = x[k][j][i + 1];
1379e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
1380e77315bcSPatrick Sanan           be = PetscPowScalar(ae, bm1);
1381e77315bcSPatrick Sanan           /* de = be * ae; */
1382e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
1383e77315bcSPatrick Sanan           ge = coef * be * (te - t0);
1384e77315bcSPatrick Sanan 
1385e77315bcSPatrick Sanan           ts = x[k][j - 1][i];
1386e77315bcSPatrick Sanan           as = 0.5 * (t0 + ts);
1387e77315bcSPatrick Sanan           bs = PetscPowScalar(as, bm1);
1388e77315bcSPatrick Sanan           /* ds = bs * as; */
1389e77315bcSPatrick Sanan           ds = PetscPowScalar(as, beta);
1390e77315bcSPatrick Sanan           gs = coef * bs * (t0 - ts);
1391e77315bcSPatrick Sanan 
1392e77315bcSPatrick Sanan           /* top down interior edge */
1393e77315bcSPatrick Sanan           if (k == 0) {
1394e77315bcSPatrick Sanan             tu = x[k + 1][j][i];
1395e77315bcSPatrick Sanan             au = 0.5 * (t0 + tu);
1396e77315bcSPatrick Sanan             bu = PetscPowScalar(au, bm1);
1397e77315bcSPatrick Sanan             /* du = bu * au; */
1398e77315bcSPatrick Sanan             du = PetscPowScalar(au, beta);
1399e77315bcSPatrick Sanan             gu = coef * bu * (tu - t0);
1400e77315bcSPatrick Sanan 
14019371c9d4SSatish Balay             c[0].k = k;
14029371c9d4SSatish Balay             c[0].j = j - 1;
14039371c9d4SSatish Balay             c[0].i = i;
1404e77315bcSPatrick Sanan             v[0]   = -hzhxdhy * (ds - gs);
14059371c9d4SSatish Balay             c[1].k = k;
14069371c9d4SSatish Balay             c[1].j = j;
14079371c9d4SSatish Balay             c[1].i = i - 1;
1408e77315bcSPatrick Sanan             v[1]   = -hyhzdhx * (dw - gw);
14099371c9d4SSatish Balay             c[2].k = k;
14109371c9d4SSatish Balay             c[2].j = j;
14119371c9d4SSatish Balay             c[2].i = i;
1412e77315bcSPatrick Sanan             v[2]   = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu);
14139371c9d4SSatish Balay             c[3].k = k;
14149371c9d4SSatish Balay             c[3].j = j;
14159371c9d4SSatish Balay             c[3].i = i + 1;
1416e77315bcSPatrick Sanan             v[3]   = -hyhzdhx * (de + ge);
14179371c9d4SSatish Balay             c[4].k = k + 1;
14189371c9d4SSatish Balay             c[4].j = j;
14199371c9d4SSatish Balay             c[4].i = i;
1420e77315bcSPatrick Sanan             v[4]   = -hxhydhz * (du + gu);
14219566063dSJacob Faibussowitsch             PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
1422e77315bcSPatrick Sanan 
1423e77315bcSPatrick Sanan           } else if (k == mz - 1) { /* top up interior edge */
1424e77315bcSPatrick Sanan 
1425e77315bcSPatrick Sanan             td = x[k - 1][j][i];
1426e77315bcSPatrick Sanan             ad = 0.5 * (t0 + td);
1427e77315bcSPatrick Sanan             bd = PetscPowScalar(ad, bm1);
1428e77315bcSPatrick Sanan             /* dd = bd * ad; */
1429e77315bcSPatrick Sanan             dd = PetscPowScalar(ad, beta);
1430e77315bcSPatrick Sanan             gd = coef * bd * (td - t0);
1431e77315bcSPatrick Sanan 
14329371c9d4SSatish Balay             c[0].k = k - 1;
14339371c9d4SSatish Balay             c[0].j = j;
14349371c9d4SSatish Balay             c[0].i = i;
1435e77315bcSPatrick Sanan             v[0]   = -hxhydhz * (dd - gd);
14369371c9d4SSatish Balay             c[1].k = k;
14379371c9d4SSatish Balay             c[1].j = j - 1;
14389371c9d4SSatish Balay             c[1].i = i;
1439e77315bcSPatrick Sanan             v[1]   = -hzhxdhy * (ds - gs);
14409371c9d4SSatish Balay             c[2].k = k;
14419371c9d4SSatish Balay             c[2].j = j;
14429371c9d4SSatish Balay             c[2].i = i - 1;
1443e77315bcSPatrick Sanan             v[2]   = -hyhzdhx * (dw - gw);
14449371c9d4SSatish Balay             c[3].k = k;
14459371c9d4SSatish Balay             c[3].j = j;
14469371c9d4SSatish Balay             c[3].i = i;
1447e77315bcSPatrick Sanan             v[3]   = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd);
14489371c9d4SSatish Balay             c[4].k = k;
14499371c9d4SSatish Balay             c[4].j = j;
14509371c9d4SSatish Balay             c[4].i = i + 1;
1451e77315bcSPatrick Sanan             v[4]   = -hyhzdhx * (de + ge);
14529566063dSJacob Faibussowitsch             PetscCall(MatSetValuesStencil(jac, 1, &row, 5, c, v, INSERT_VALUES));
1453e77315bcSPatrick Sanan 
1454e77315bcSPatrick Sanan           } else { /* top interior plane */
1455e77315bcSPatrick Sanan 
1456e77315bcSPatrick Sanan             tu = x[k + 1][j][i];
1457e77315bcSPatrick Sanan             au = 0.5 * (t0 + tu);
1458e77315bcSPatrick Sanan             bu = PetscPowScalar(au, bm1);
1459e77315bcSPatrick Sanan             /* du = bu * au; */
1460e77315bcSPatrick Sanan             du = PetscPowScalar(au, beta);
1461e77315bcSPatrick Sanan             gu = coef * bu * (tu - t0);
1462e77315bcSPatrick Sanan 
1463e77315bcSPatrick Sanan             td = x[k - 1][j][i];
1464e77315bcSPatrick Sanan             ad = 0.5 * (t0 + td);
1465e77315bcSPatrick Sanan             bd = PetscPowScalar(ad, bm1);
1466e77315bcSPatrick Sanan             /* dd = bd * ad; */
1467e77315bcSPatrick Sanan             dd = PetscPowScalar(ad, beta);
1468e77315bcSPatrick Sanan             gd = coef * bd * (td - t0);
1469e77315bcSPatrick Sanan 
14709371c9d4SSatish Balay             c[0].k = k - 1;
14719371c9d4SSatish Balay             c[0].j = j;
14729371c9d4SSatish Balay             c[0].i = i;
1473e77315bcSPatrick Sanan             v[0]   = -hxhydhz * (dd - gd);
14749371c9d4SSatish Balay             c[1].k = k;
14759371c9d4SSatish Balay             c[1].j = j - 1;
14769371c9d4SSatish Balay             c[1].i = i;
1477e77315bcSPatrick Sanan             v[1]   = -hzhxdhy * (ds - gs);
14789371c9d4SSatish Balay             c[2].k = k;
14799371c9d4SSatish Balay             c[2].j = j;
14809371c9d4SSatish Balay             c[2].i = i - 1;
1481e77315bcSPatrick Sanan             v[2]   = -hyhzdhx * (dw - gw);
14829371c9d4SSatish Balay             c[3].k = k;
14839371c9d4SSatish Balay             c[3].j = j;
14849371c9d4SSatish Balay             c[3].i = i;
1485e77315bcSPatrick Sanan             v[3]   = hzhxdhy * (ds + gs) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + du + gd - gu);
14869371c9d4SSatish Balay             c[4].k = k;
14879371c9d4SSatish Balay             c[4].j = j;
14889371c9d4SSatish Balay             c[4].i = i + 1;
1489e77315bcSPatrick Sanan             v[4]   = -hyhzdhx * (de + ge);
14909371c9d4SSatish Balay             c[5].k = k + 1;
14919371c9d4SSatish Balay             c[5].j = j;
14929371c9d4SSatish Balay             c[5].i = i;
1493e77315bcSPatrick Sanan             v[5]   = -hxhydhz * (du + gu);
14949566063dSJacob Faibussowitsch             PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES));
1495e77315bcSPatrick Sanan           }
1496e77315bcSPatrick Sanan 
1497e77315bcSPatrick Sanan         } else if (k == 0) {
1498e77315bcSPatrick Sanan           /* down interior plane */
1499e77315bcSPatrick Sanan 
1500e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
1501e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
1502e77315bcSPatrick Sanan           bw = PetscPowScalar(aw, bm1);
1503e77315bcSPatrick Sanan           /* dw = bw * aw */
1504e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
1505e77315bcSPatrick Sanan           gw = coef * bw * (t0 - tw);
1506e77315bcSPatrick Sanan 
1507e77315bcSPatrick Sanan           te = x[k][j][i + 1];
1508e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
1509e77315bcSPatrick Sanan           be = PetscPowScalar(ae, bm1);
1510e77315bcSPatrick Sanan           /* de = be * ae; */
1511e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
1512e77315bcSPatrick Sanan           ge = coef * be * (te - t0);
1513e77315bcSPatrick Sanan 
1514e77315bcSPatrick Sanan           ts = x[k][j - 1][i];
1515e77315bcSPatrick Sanan           as = 0.5 * (t0 + ts);
1516e77315bcSPatrick Sanan           bs = PetscPowScalar(as, bm1);
1517e77315bcSPatrick Sanan           /* ds = bs * as; */
1518e77315bcSPatrick Sanan           ds = PetscPowScalar(as, beta);
1519e77315bcSPatrick Sanan           gs = coef * bs * (t0 - ts);
1520e77315bcSPatrick Sanan 
1521e77315bcSPatrick Sanan           tn = x[k][j + 1][i];
1522e77315bcSPatrick Sanan           an = 0.5 * (t0 + tn);
1523e77315bcSPatrick Sanan           bn = PetscPowScalar(an, bm1);
1524e77315bcSPatrick Sanan           /* dn = bn * an; */
1525e77315bcSPatrick Sanan           dn = PetscPowScalar(an, beta);
1526e77315bcSPatrick Sanan           gn = coef * bn * (tn - t0);
1527e77315bcSPatrick Sanan 
1528e77315bcSPatrick Sanan           tu = x[k + 1][j][i];
1529e77315bcSPatrick Sanan           au = 0.5 * (t0 + tu);
1530e77315bcSPatrick Sanan           bu = PetscPowScalar(au, bm1);
1531e77315bcSPatrick Sanan           /* du = bu * au; */
1532e77315bcSPatrick Sanan           du = PetscPowScalar(au, beta);
1533e77315bcSPatrick Sanan           gu = coef * bu * (tu - t0);
1534e77315bcSPatrick Sanan 
15359371c9d4SSatish Balay           c[0].k = k;
15369371c9d4SSatish Balay           c[0].j = j - 1;
15379371c9d4SSatish Balay           c[0].i = i;
1538e77315bcSPatrick Sanan           v[0]   = -hzhxdhy * (ds - gs);
15399371c9d4SSatish Balay           c[1].k = k;
15409371c9d4SSatish Balay           c[1].j = j;
15419371c9d4SSatish Balay           c[1].i = i - 1;
1542e77315bcSPatrick Sanan           v[1]   = -hyhzdhx * (dw - gw);
15439371c9d4SSatish Balay           c[2].k = k;
15449371c9d4SSatish Balay           c[2].j = j;
15459371c9d4SSatish Balay           c[2].i = i;
1546e77315bcSPatrick Sanan           v[2]   = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (du - gu);
15479371c9d4SSatish Balay           c[3].k = k;
15489371c9d4SSatish Balay           c[3].j = j;
15499371c9d4SSatish Balay           c[3].i = i + 1;
1550e77315bcSPatrick Sanan           v[3]   = -hyhzdhx * (de + ge);
15519371c9d4SSatish Balay           c[4].k = k;
15529371c9d4SSatish Balay           c[4].j = j + 1;
15539371c9d4SSatish Balay           c[4].i = i;
1554e77315bcSPatrick Sanan           v[4]   = -hzhxdhy * (dn + gn);
15559371c9d4SSatish Balay           c[5].k = k + 1;
15569371c9d4SSatish Balay           c[5].j = j;
15579371c9d4SSatish Balay           c[5].i = i;
1558e77315bcSPatrick Sanan           v[5]   = -hxhydhz * (du + gu);
15599566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES));
1560e77315bcSPatrick Sanan 
1561e77315bcSPatrick Sanan         } else if (k == mz - 1) {
1562e77315bcSPatrick Sanan           /* up interior plane */
1563e77315bcSPatrick Sanan 
1564e77315bcSPatrick Sanan           tw = x[k][j][i - 1];
1565e77315bcSPatrick Sanan           aw = 0.5 * (t0 + tw);
1566e77315bcSPatrick Sanan           bw = PetscPowScalar(aw, bm1);
1567e77315bcSPatrick Sanan           /* dw = bw * aw */
1568e77315bcSPatrick Sanan           dw = PetscPowScalar(aw, beta);
1569e77315bcSPatrick Sanan           gw = coef * bw * (t0 - tw);
1570e77315bcSPatrick Sanan 
1571e77315bcSPatrick Sanan           te = x[k][j][i + 1];
1572e77315bcSPatrick Sanan           ae = 0.5 * (t0 + te);
1573e77315bcSPatrick Sanan           be = PetscPowScalar(ae, bm1);
1574e77315bcSPatrick Sanan           /* de = be * ae; */
1575e77315bcSPatrick Sanan           de = PetscPowScalar(ae, beta);
1576e77315bcSPatrick Sanan           ge = coef * be * (te - t0);
1577e77315bcSPatrick Sanan 
1578e77315bcSPatrick Sanan           ts = x[k][j - 1][i];
1579e77315bcSPatrick Sanan           as = 0.5 * (t0 + ts);
1580e77315bcSPatrick Sanan           bs = PetscPowScalar(as, bm1);
1581e77315bcSPatrick Sanan           /* ds = bs * as; */
1582e77315bcSPatrick Sanan           ds = PetscPowScalar(as, beta);
1583e77315bcSPatrick Sanan           gs = coef * bs * (t0 - ts);
1584e77315bcSPatrick Sanan 
1585e77315bcSPatrick Sanan           tn = x[k][j + 1][i];
1586e77315bcSPatrick Sanan           an = 0.5 * (t0 + tn);
1587e77315bcSPatrick Sanan           bn = PetscPowScalar(an, bm1);
1588e77315bcSPatrick Sanan           /* dn = bn * an; */
1589e77315bcSPatrick Sanan           dn = PetscPowScalar(an, beta);
1590e77315bcSPatrick Sanan           gn = coef * bn * (tn - t0);
1591e77315bcSPatrick Sanan 
1592e77315bcSPatrick Sanan           td = x[k - 1][j][i];
1593e77315bcSPatrick Sanan           ad = 0.5 * (t0 + td);
1594e77315bcSPatrick Sanan           bd = PetscPowScalar(ad, bm1);
1595e77315bcSPatrick Sanan           /* dd = bd * ad; */
1596e77315bcSPatrick Sanan           dd = PetscPowScalar(ad, beta);
1597e77315bcSPatrick Sanan           gd = coef * bd * (t0 - td);
1598e77315bcSPatrick Sanan 
15999371c9d4SSatish Balay           c[0].k = k - 1;
16009371c9d4SSatish Balay           c[0].j = j;
16019371c9d4SSatish Balay           c[0].i = i;
1602e77315bcSPatrick Sanan           v[0]   = -hxhydhz * (dd - gd);
16039371c9d4SSatish Balay           c[1].k = k;
16049371c9d4SSatish Balay           c[1].j = j - 1;
16059371c9d4SSatish Balay           c[1].i = i;
1606e77315bcSPatrick Sanan           v[1]   = -hzhxdhy * (ds - gs);
16079371c9d4SSatish Balay           c[2].k = k;
16089371c9d4SSatish Balay           c[2].j = j;
16099371c9d4SSatish Balay           c[2].i = i - 1;
1610e77315bcSPatrick Sanan           v[2]   = -hyhzdhx * (dw - gw);
16119371c9d4SSatish Balay           c[3].k = k;
16129371c9d4SSatish Balay           c[3].j = j;
16139371c9d4SSatish Balay           c[3].i = i;
1614e77315bcSPatrick Sanan           v[3]   = hzhxdhy * (ds + dn + gs - gn) + hyhzdhx * (dw + de + gw - ge) + hxhydhz * (dd + gd);
16159371c9d4SSatish Balay           c[4].k = k;
16169371c9d4SSatish Balay           c[4].j = j;
16179371c9d4SSatish Balay           c[4].i = i + 1;
1618e77315bcSPatrick Sanan           v[4]   = -hyhzdhx * (de + ge);
16199371c9d4SSatish Balay           c[5].k = k;
16209371c9d4SSatish Balay           c[5].j = j + 1;
16219371c9d4SSatish Balay           c[5].i = i;
1622e77315bcSPatrick Sanan           v[5]   = -hzhxdhy * (dn + gn);
16239566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(jac, 1, &row, 6, c, v, INSERT_VALUES));
1624e77315bcSPatrick Sanan         }
1625e77315bcSPatrick Sanan       }
1626e77315bcSPatrick Sanan     }
1627e77315bcSPatrick Sanan   }
16289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(jac, MAT_FINAL_ASSEMBLY));
16299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(jac, MAT_FINAL_ASSEMBLY));
1630e77315bcSPatrick Sanan   if (jac != J) {
16319566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
16329566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
1633e77315bcSPatrick Sanan   }
16349566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, localX, &x));
16359566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
1636e77315bcSPatrick Sanan 
16379566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops((41.0 + 8.0 * POWFLOP) * xm * ym));
16383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1639e77315bcSPatrick Sanan }
1640e77315bcSPatrick Sanan 
TestConvergence(SNES snes,PETSC_UNUSED PetscInt it,PETSC_UNUSED PetscReal xnorm,PETSC_UNUSED PetscReal snorm,PETSC_UNUSED PetscReal fnorm,SNESConvergedReason * reason,PetscCtx ctx)1641*2a8381b2SBarry Smith PetscErrorCode TestConvergence(SNES snes, PETSC_UNUSED PetscInt it, PETSC_UNUSED PetscReal xnorm, PETSC_UNUSED PetscReal snorm, PETSC_UNUSED PetscReal fnorm, SNESConvergedReason *reason, PetscCtx ctx)
1642379ef8d2SAlexander {
1643379ef8d2SAlexander   AppCtx *user = (AppCtx *)ctx;
1644379ef8d2SAlexander 
1645379ef8d2SAlexander   PetscFunctionBeginUser;
1646379ef8d2SAlexander   if (user->converged) *reason = SNES_CONVERGED_USER;
1647379ef8d2SAlexander   else *reason = SNES_DIVERGED_USER;
1648379ef8d2SAlexander   PetscFunctionReturn(PETSC_SUCCESS);
1649379ef8d2SAlexander }
1650379ef8d2SAlexander 
1651e77315bcSPatrick Sanan /*TEST
1652e77315bcSPatrick Sanan 
1653e77315bcSPatrick Sanan    test:
1654e77315bcSPatrick Sanan       nsize: 4
1655e77315bcSPatrick Sanan       args: -snes_monitor_short -pc_mg_type full -ksp_type fgmres -pc_type mg -snes_view -pc_mg_levels 2 -pc_mg_galerkin pmat
1656e77315bcSPatrick Sanan       requires: !single
1657e77315bcSPatrick Sanan 
1658379ef8d2SAlexander    test:
1659379ef8d2SAlexander       suffix: 2
1660379ef8d2SAlexander       args: -snes_converged_reason -use_convergence_test -mark_converged
1661379ef8d2SAlexander 
1662379ef8d2SAlexander    test:
1663379ef8d2SAlexander       suffix: 3
1664379ef8d2SAlexander       args: -snes_converged_reason -use_convergence_test -mark_converged 0
1665379ef8d2SAlexander 
1666e77315bcSPatrick Sanan TEST*/
1667