xref: /petsc/src/snes/tests/ex20.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1e77315bcSPatrick Sanan 
2e77315bcSPatrick Sanan static char help[] ="Nonlinear Radiative Transport PDE with multigrid in 3d.\n\
3e77315bcSPatrick Sanan Uses 3-dimensional distributed arrays.\n\
4e77315bcSPatrick Sanan A 3-dim simplified Radiative Transport test problem is used, with analytic Jacobian. \n\
5e77315bcSPatrick Sanan \n\
6e77315bcSPatrick Sanan   Solves the linear systems via multilevel methods \n\
7e77315bcSPatrick Sanan \n\
8e77315bcSPatrick Sanan The command line\n\
9e77315bcSPatrick Sanan options are:\n\
10e77315bcSPatrick Sanan   -tleft <tl>, where <tl> indicates the left Diriclet BC \n\
11e77315bcSPatrick Sanan   -tright <tr>, where <tr> indicates the right Diriclet BC \n\
12e77315bcSPatrick Sanan   -beta <beta>, where <beta> indicates the exponent in T \n\n";
13e77315bcSPatrick Sanan 
14e77315bcSPatrick Sanan /*T
15e77315bcSPatrick Sanan    Concepts: SNES^solving a system of nonlinear equations
16e77315bcSPatrick Sanan    Concepts: DM^using distributed arrays
17e77315bcSPatrick Sanan    Concepts: multigrid;
18e77315bcSPatrick Sanan    Processors: n
19e77315bcSPatrick Sanan T*/
20e77315bcSPatrick Sanan 
21e77315bcSPatrick Sanan /*
22e77315bcSPatrick Sanan 
23e77315bcSPatrick Sanan     This example models the partial differential equation
24e77315bcSPatrick Sanan 
25e77315bcSPatrick Sanan          - Div(alpha* T^beta (GRAD T)) = 0.
26e77315bcSPatrick Sanan 
27e77315bcSPatrick Sanan     where beta = 2.5 and alpha = 1.0
28e77315bcSPatrick Sanan 
29e77315bcSPatrick Sanan     BC: T_left = 1.0, T_right = 0.1, dT/dn_top = dTdn_bottom = dT/dn_up = dT/dn_down = 0.
30e77315bcSPatrick Sanan 
31e77315bcSPatrick Sanan     in the unit square, which is uniformly discretized in each of x and
32e77315bcSPatrick Sanan     y in this simple encoding.  The degrees of freedom are cell centered.
33e77315bcSPatrick Sanan 
34e77315bcSPatrick Sanan     A finite volume approximation with the usual 7-point stencil
35e77315bcSPatrick Sanan     is used to discretize the boundary value problem to obtain a
36e77315bcSPatrick Sanan     nonlinear system of equations.
37e77315bcSPatrick Sanan 
38e77315bcSPatrick Sanan     This code was contributed by Nickolas Jovanovic based on ex18.c
39e77315bcSPatrick Sanan 
40e77315bcSPatrick Sanan */
41e77315bcSPatrick Sanan 
42e77315bcSPatrick Sanan #include <petscsnes.h>
43e77315bcSPatrick Sanan #include <petscdm.h>
44e77315bcSPatrick Sanan #include <petscdmda.h>
45e77315bcSPatrick Sanan 
46e77315bcSPatrick Sanan /* User-defined application context */
47e77315bcSPatrick Sanan 
48e77315bcSPatrick Sanan typedef struct {
49e77315bcSPatrick Sanan   PetscReal tleft,tright;     /* Dirichlet boundary conditions */
50e77315bcSPatrick Sanan   PetscReal beta,bm1,coef;    /* nonlinear diffusivity parameterizations */
51e77315bcSPatrick Sanan } AppCtx;
52e77315bcSPatrick Sanan 
53e77315bcSPatrick Sanan #define POWFLOP 5 /* assume a pow() takes five flops */
54e77315bcSPatrick Sanan 
55e77315bcSPatrick Sanan extern PetscErrorCode FormInitialGuess(SNES,Vec,void*);
56e77315bcSPatrick Sanan extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
57e77315bcSPatrick Sanan extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
58e77315bcSPatrick Sanan 
59e77315bcSPatrick Sanan int main(int argc,char **argv)
60e77315bcSPatrick Sanan {
61e77315bcSPatrick Sanan   SNES           snes;
62e77315bcSPatrick Sanan   AppCtx         user;
63e77315bcSPatrick Sanan   PetscInt       its,lits;
64e77315bcSPatrick Sanan   PetscReal      litspit;
65e77315bcSPatrick Sanan   DM             da;
66e77315bcSPatrick Sanan 
67*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
68e77315bcSPatrick Sanan   /* set problem parameters */
69e77315bcSPatrick Sanan   user.tleft  = 1.0;
70e77315bcSPatrick Sanan   user.tright = 0.1;
71e77315bcSPatrick Sanan   user.beta   = 2.5;
725f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tleft",&user.tleft,NULL));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tright",&user.tright,NULL));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-beta",&user.beta,NULL));
75e77315bcSPatrick Sanan   user.bm1    = user.beta - 1.0;
76e77315bcSPatrick Sanan   user.coef   = user.beta/2.0;
77e77315bcSPatrick Sanan 
78e77315bcSPatrick Sanan   /*
79e77315bcSPatrick Sanan       Set the DMDA (grid structure) for the grids.
80e77315bcSPatrick Sanan   */
815f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(da,&user));
85e77315bcSPatrick Sanan 
86e77315bcSPatrick Sanan   /*
87e77315bcSPatrick Sanan      Create the nonlinear solver
88e77315bcSPatrick Sanan   */
895f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes,da));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,NULL,FormFunction,&user));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,NULL,NULL,FormJacobian,&user));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetComputeInitialGuess(snes,FormInitialGuess,NULL));
95e77315bcSPatrick Sanan 
965f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,NULL));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetLinearSolveIterations(snes,&lits));
99e77315bcSPatrick Sanan   litspit = ((PetscReal)lits)/((PetscReal)its);
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of Linear iterations = %D\n",lits));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Average Linear its / SNES = %e\n",(double)litspit));
103e77315bcSPatrick Sanan 
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
106*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
107*b122ec5aSJacob Faibussowitsch   return 0;
108e77315bcSPatrick Sanan }
109e77315bcSPatrick Sanan /* --------------------  Form initial approximation ----------------- */
110e77315bcSPatrick Sanan PetscErrorCode FormInitialGuess(SNES snes,Vec X,void *ctx)
111e77315bcSPatrick Sanan {
112e77315bcSPatrick Sanan   AppCtx         *user;
113e77315bcSPatrick Sanan   PetscInt       i,j,k,xs,ys,xm,ym,zs,zm;
114e77315bcSPatrick Sanan   PetscScalar    ***x;
115e77315bcSPatrick Sanan   DM             da;
116e77315bcSPatrick Sanan 
117e77315bcSPatrick Sanan   PetscFunctionBeginUser;
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(da,&user));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
122e77315bcSPatrick Sanan 
123e77315bcSPatrick Sanan   /* Compute initial guess */
124e77315bcSPatrick Sanan   for (k=zs; k<zs+zm; k++) {
125e77315bcSPatrick Sanan     for (j=ys; j<ys+ym; j++) {
126e77315bcSPatrick Sanan       for (i=xs; i<xs+xm; i++) {
127e77315bcSPatrick Sanan         x[k][j][i] = user->tleft;
128e77315bcSPatrick Sanan       }
129e77315bcSPatrick Sanan     }
130e77315bcSPatrick Sanan   }
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
132e77315bcSPatrick Sanan   PetscFunctionReturn(0);
133e77315bcSPatrick Sanan }
134e77315bcSPatrick Sanan /* --------------------  Evaluate Function F(x) --------------------- */
135e77315bcSPatrick Sanan PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
136e77315bcSPatrick Sanan {
137e77315bcSPatrick Sanan   AppCtx         *user = (AppCtx*)ptr;
138e77315bcSPatrick Sanan   PetscInt       i,j,k,mx,my,mz,xs,ys,zs,xm,ym,zm;
139e77315bcSPatrick Sanan   PetscScalar    zero = 0.0,one = 1.0;
140e77315bcSPatrick Sanan   PetscScalar    hx,hy,hz,hxhydhz,hyhzdhx,hzhxdhy;
141e77315bcSPatrick 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;
142e77315bcSPatrick Sanan   PetscScalar    tleft,tright,beta,td,ad,dd,fd=0.0,tu,au,du=0.0,fu=0.0;
143e77315bcSPatrick Sanan   PetscScalar    ***x,***f;
144e77315bcSPatrick Sanan   Vec            localX;
145e77315bcSPatrick Sanan   DM             da;
146e77315bcSPatrick Sanan 
147e77315bcSPatrick Sanan   PetscFunctionBeginUser;
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localX));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
151e77315bcSPatrick Sanan   hx      = one/(PetscReal)(mx-1);  hy    = one/(PetscReal)(my-1);  hz = one/(PetscReal)(mz-1);
152e77315bcSPatrick Sanan   hxhydhz = hx*hy/hz;   hyhzdhx = hy*hz/hx;   hzhxdhy = hz*hx/hy;
153e77315bcSPatrick Sanan   tleft   = user->tleft;         tright = user->tright;
154e77315bcSPatrick Sanan   beta    = user->beta;
155e77315bcSPatrick Sanan 
156e77315bcSPatrick Sanan   /* Get ghost points */
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,localX,&x));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,F,&f));
162e77315bcSPatrick Sanan 
163e77315bcSPatrick Sanan   /* Evaluate function */
164e77315bcSPatrick Sanan   for (k=zs; k<zs+zm; k++) {
165e77315bcSPatrick Sanan     for (j=ys; j<ys+ym; j++) {
166e77315bcSPatrick Sanan       for (i=xs; i<xs+xm; i++) {
167e77315bcSPatrick Sanan         t0 = x[k][j][i];
168e77315bcSPatrick Sanan 
169e77315bcSPatrick Sanan         if (i > 0 && i < mx-1 && j > 0 && j < my-1 && k > 0 && k < mz-1) {
170e77315bcSPatrick Sanan 
171e77315bcSPatrick Sanan           /* general interior volume */
172e77315bcSPatrick Sanan 
173e77315bcSPatrick Sanan           tw = x[k][j][i-1];
174e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
175e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
176e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
177e77315bcSPatrick Sanan 
178e77315bcSPatrick Sanan           te = x[k][j][i+1];
179e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
180e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
181e77315bcSPatrick Sanan           fe = de*(te - t0);
182e77315bcSPatrick Sanan 
183e77315bcSPatrick Sanan           ts = x[k][j-1][i];
184e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
185e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
186e77315bcSPatrick Sanan           fs = ds*(t0 - ts);
187e77315bcSPatrick Sanan 
188e77315bcSPatrick Sanan           tn = x[k][j+1][i];
189e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
190e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
191e77315bcSPatrick Sanan           fn = dn*(tn - t0);
192e77315bcSPatrick Sanan 
193e77315bcSPatrick Sanan           td = x[k-1][j][i];
194e77315bcSPatrick Sanan           ad = 0.5*(t0 + td);
195e77315bcSPatrick Sanan           dd = PetscPowScalar(ad,beta);
196e77315bcSPatrick Sanan           fd = dd*(t0 - td);
197e77315bcSPatrick Sanan 
198e77315bcSPatrick Sanan           tu = x[k+1][j][i];
199e77315bcSPatrick Sanan           au = 0.5*(t0 + tu);
200e77315bcSPatrick Sanan           du = PetscPowScalar(au,beta);
201e77315bcSPatrick Sanan           fu = du*(tu - t0);
202e77315bcSPatrick Sanan 
203e77315bcSPatrick Sanan         } else if (i == 0) {
204e77315bcSPatrick Sanan 
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 
254e77315bcSPatrick Sanan           /* right-hand (east) boundary */
255e77315bcSPatrick Sanan           tw = x[k][j][i-1];
256e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
257e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
258e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
259e77315bcSPatrick Sanan 
260e77315bcSPatrick Sanan           te = tright;
261e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
262e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
263e77315bcSPatrick Sanan           fe = de*(te - t0);
264e77315bcSPatrick Sanan 
265e77315bcSPatrick Sanan           if (j > 0) {
266e77315bcSPatrick Sanan             ts = x[k][j-1][i];
267e77315bcSPatrick Sanan             as = 0.5*(t0 + ts);
268e77315bcSPatrick Sanan             ds = PetscPowScalar(as,beta);
269e77315bcSPatrick Sanan             fs = ds*(t0 - ts);
270e77315bcSPatrick Sanan           } else {
271e77315bcSPatrick Sanan             fs = zero;
272e77315bcSPatrick Sanan           }
273e77315bcSPatrick Sanan 
274e77315bcSPatrick Sanan           if (j < my-1) {
275e77315bcSPatrick Sanan             tn = x[k][j+1][i];
276e77315bcSPatrick Sanan             an = 0.5*(t0 + tn);
277e77315bcSPatrick Sanan             dn = PetscPowScalar(an,beta);
278e77315bcSPatrick Sanan             fn = dn*(tn - t0);
279e77315bcSPatrick Sanan           } else {
280e77315bcSPatrick Sanan             fn = zero;
281e77315bcSPatrick Sanan           }
282e77315bcSPatrick Sanan 
283e77315bcSPatrick Sanan           if (k > 0) {
284e77315bcSPatrick Sanan             td = x[k-1][j][i];
285e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
286e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
287e77315bcSPatrick Sanan             fd = dd*(t0 - td);
288e77315bcSPatrick Sanan           } else {
289e77315bcSPatrick Sanan             fd = zero;
290e77315bcSPatrick Sanan           }
291e77315bcSPatrick Sanan 
292e77315bcSPatrick Sanan           if (k < mz-1) {
293e77315bcSPatrick Sanan             tu = x[k+1][j][i];
294e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
295e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
296e77315bcSPatrick Sanan             fu = du*(tu - t0);
297e77315bcSPatrick Sanan           } else {
298e77315bcSPatrick Sanan             fu = zero;
299e77315bcSPatrick Sanan           }
300e77315bcSPatrick Sanan 
301e77315bcSPatrick Sanan         } else if (j == 0) {
302e77315bcSPatrick Sanan 
303e77315bcSPatrick Sanan           /* bottom (south) boundary, and i <> 0 or mx-1 */
304e77315bcSPatrick Sanan           tw = x[k][j][i-1];
305e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
306e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
307e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
308e77315bcSPatrick Sanan 
309e77315bcSPatrick Sanan           te = x[k][j][i+1];
310e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
311e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
312e77315bcSPatrick Sanan           fe = de*(te - t0);
313e77315bcSPatrick Sanan 
314e77315bcSPatrick Sanan           fs = zero;
315e77315bcSPatrick Sanan 
316e77315bcSPatrick Sanan           tn = x[k][j+1][i];
317e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
318e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
319e77315bcSPatrick Sanan           fn = dn*(tn - t0);
320e77315bcSPatrick Sanan 
321e77315bcSPatrick Sanan           if (k > 0) {
322e77315bcSPatrick Sanan             td = x[k-1][j][i];
323e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
324e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
325e77315bcSPatrick Sanan             fd = dd*(t0 - td);
326e77315bcSPatrick Sanan           } else {
327e77315bcSPatrick Sanan             fd = zero;
328e77315bcSPatrick Sanan           }
329e77315bcSPatrick Sanan 
330e77315bcSPatrick Sanan           if (k < mz-1) {
331e77315bcSPatrick Sanan             tu = x[k+1][j][i];
332e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
333e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
334e77315bcSPatrick Sanan             fu = du*(tu - t0);
335e77315bcSPatrick Sanan           } else {
336e77315bcSPatrick Sanan             fu = zero;
337e77315bcSPatrick Sanan           }
338e77315bcSPatrick Sanan 
339e77315bcSPatrick Sanan         } else if (j == my-1) {
340e77315bcSPatrick Sanan 
341e77315bcSPatrick Sanan           /* top (north) boundary, and i <> 0 or mx-1 */
342e77315bcSPatrick Sanan           tw = x[k][j][i-1];
343e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
344e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
345e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
346e77315bcSPatrick Sanan 
347e77315bcSPatrick Sanan           te = x[k][j][i+1];
348e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
349e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
350e77315bcSPatrick Sanan           fe = de*(te - t0);
351e77315bcSPatrick Sanan 
352e77315bcSPatrick Sanan           ts = x[k][j-1][i];
353e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
354e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
355e77315bcSPatrick Sanan           fs = ds*(t0 - ts);
356e77315bcSPatrick Sanan 
357e77315bcSPatrick Sanan           fn = zero;
358e77315bcSPatrick Sanan 
359e77315bcSPatrick Sanan           if (k > 0) {
360e77315bcSPatrick Sanan             td = x[k-1][j][i];
361e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
362e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
363e77315bcSPatrick Sanan             fd = dd*(t0 - td);
364e77315bcSPatrick Sanan           } else {
365e77315bcSPatrick Sanan             fd = zero;
366e77315bcSPatrick Sanan           }
367e77315bcSPatrick Sanan 
368e77315bcSPatrick Sanan           if (k < mz-1) {
369e77315bcSPatrick Sanan             tu = x[k+1][j][i];
370e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
371e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
372e77315bcSPatrick Sanan             fu = du*(tu - t0);
373e77315bcSPatrick Sanan           } else {
374e77315bcSPatrick Sanan             fu = zero;
375e77315bcSPatrick Sanan           }
376e77315bcSPatrick Sanan 
377e77315bcSPatrick Sanan         } else if (k == 0) {
378e77315bcSPatrick Sanan 
379e77315bcSPatrick Sanan           /* down boundary (interior only) */
380e77315bcSPatrick Sanan           tw = x[k][j][i-1];
381e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
382e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
383e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
384e77315bcSPatrick Sanan 
385e77315bcSPatrick Sanan           te = x[k][j][i+1];
386e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
387e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
388e77315bcSPatrick Sanan           fe = de*(te - t0);
389e77315bcSPatrick Sanan 
390e77315bcSPatrick Sanan           ts = x[k][j-1][i];
391e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
392e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
393e77315bcSPatrick Sanan           fs = ds*(t0 - ts);
394e77315bcSPatrick Sanan 
395e77315bcSPatrick Sanan           tn = x[k][j+1][i];
396e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
397e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
398e77315bcSPatrick Sanan           fn = dn*(tn - t0);
399e77315bcSPatrick Sanan 
400e77315bcSPatrick Sanan           fd = zero;
401e77315bcSPatrick Sanan 
402e77315bcSPatrick Sanan           tu = x[k+1][j][i];
403e77315bcSPatrick Sanan           au = 0.5*(t0 + tu);
404e77315bcSPatrick Sanan           du = PetscPowScalar(au,beta);
405e77315bcSPatrick Sanan           fu = du*(tu - t0);
406e77315bcSPatrick Sanan 
407e77315bcSPatrick Sanan         } else if (k == mz-1) {
408e77315bcSPatrick Sanan 
409e77315bcSPatrick Sanan           /* up boundary (interior only) */
410e77315bcSPatrick Sanan           tw = x[k][j][i-1];
411e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
412e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
413e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
414e77315bcSPatrick Sanan 
415e77315bcSPatrick Sanan           te = x[k][j][i+1];
416e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
417e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
418e77315bcSPatrick Sanan           fe = de*(te - t0);
419e77315bcSPatrick Sanan 
420e77315bcSPatrick Sanan           ts = x[k][j-1][i];
421e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
422e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
423e77315bcSPatrick Sanan           fs = ds*(t0 - ts);
424e77315bcSPatrick Sanan 
425e77315bcSPatrick Sanan           tn = x[k][j+1][i];
426e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
427e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
428e77315bcSPatrick Sanan           fn = dn*(tn - t0);
429e77315bcSPatrick Sanan 
430e77315bcSPatrick Sanan           td = x[k-1][j][i];
431e77315bcSPatrick Sanan           ad = 0.5*(t0 + td);
432e77315bcSPatrick Sanan           dd = PetscPowScalar(ad,beta);
433e77315bcSPatrick Sanan           fd = dd*(t0 - td);
434e77315bcSPatrick Sanan 
435e77315bcSPatrick Sanan           fu = zero;
436e77315bcSPatrick Sanan         }
437e77315bcSPatrick Sanan 
438e77315bcSPatrick Sanan         f[k][j][i] = -hyhzdhx*(fe-fw) - hzhxdhy*(fn-fs) - hxhydhz*(fu-fd);
439e77315bcSPatrick Sanan       }
440e77315bcSPatrick Sanan     }
441e77315bcSPatrick Sanan   }
4425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,localX,&x));
4435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,F,&f));
4445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localX));
4455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops((22.0 + 4.0*POWFLOP)*ym*xm));
446e77315bcSPatrick Sanan   PetscFunctionReturn(0);
447e77315bcSPatrick Sanan }
448e77315bcSPatrick Sanan /* --------------------  Evaluate Jacobian F(x) --------------------- */
449e77315bcSPatrick Sanan PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr)
450e77315bcSPatrick Sanan {
451e77315bcSPatrick Sanan   AppCtx         *user = (AppCtx*)ptr;
452e77315bcSPatrick Sanan   PetscInt       i,j,k,mx,my,mz,xs,ys,zs,xm,ym,zm;
453e77315bcSPatrick Sanan   PetscScalar    one = 1.0;
454e77315bcSPatrick Sanan   PetscScalar    hx,hy,hz,hxhydhz,hyhzdhx,hzhxdhy;
455e77315bcSPatrick Sanan   PetscScalar    t0,tn,ts,te,tw,an,as,ae,aw,dn,ds,de,dw;
456e77315bcSPatrick Sanan   PetscScalar    tleft,tright,beta,td,ad,dd,tu,au,du,v[7],bm1,coef;
457e77315bcSPatrick Sanan   PetscScalar    ***x,bn,bs,be,bw,bu,bd,gn,gs,ge,gw,gu,gd;
458e77315bcSPatrick Sanan   Vec            localX;
459e77315bcSPatrick Sanan   MatStencil     c[7],row;
460e77315bcSPatrick Sanan   DM             da;
461e77315bcSPatrick Sanan 
462e77315bcSPatrick Sanan   PetscFunctionBeginUser;
4635f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
4645f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localX));
4655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
466e77315bcSPatrick Sanan   hx      = one/(PetscReal)(mx-1);  hy    = one/(PetscReal)(my-1);  hz = one/(PetscReal)(mz-1);
467e77315bcSPatrick Sanan   hxhydhz = hx*hy/hz;   hyhzdhx = hy*hz/hx;   hzhxdhy = hz*hx/hy;
468e77315bcSPatrick Sanan   tleft   = user->tleft;         tright = user->tright;
469e77315bcSPatrick Sanan   beta    = user->beta;          bm1    = user->bm1;              coef = user->coef;
470e77315bcSPatrick Sanan 
471e77315bcSPatrick Sanan   /* Get ghost points */
4725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
4735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
4745f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
4755f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,localX,&x));
476e77315bcSPatrick Sanan 
477e77315bcSPatrick Sanan   /* Evaluate Jacobian of function */
478e77315bcSPatrick Sanan   for (k=zs; k<zs+zm; k++) {
479e77315bcSPatrick Sanan     for (j=ys; j<ys+ym; j++) {
480e77315bcSPatrick Sanan       for (i=xs; i<xs+xm; i++) {
481e77315bcSPatrick Sanan         t0    = x[k][j][i];
482e77315bcSPatrick Sanan         row.k = k; row.j = j; row.i = i;
483e77315bcSPatrick Sanan         if (i > 0 && i < mx-1 && j > 0 && j < my-1 && k > 0 && k < mz-1) {
484e77315bcSPatrick Sanan 
485e77315bcSPatrick Sanan           /* general interior volume */
486e77315bcSPatrick Sanan 
487e77315bcSPatrick Sanan           tw = x[k][j][i-1];
488e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
489e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
490e77315bcSPatrick Sanan           /* dw = bw * aw */
491e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
492e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
493e77315bcSPatrick Sanan 
494e77315bcSPatrick Sanan           te = x[k][j][i+1];
495e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
496e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
497e77315bcSPatrick Sanan           /* de = be * ae; */
498e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
499e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
500e77315bcSPatrick Sanan 
501e77315bcSPatrick Sanan           ts = x[k][j-1][i];
502e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
503e77315bcSPatrick Sanan           bs = PetscPowScalar(as,bm1);
504e77315bcSPatrick Sanan           /* ds = bs * as; */
505e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
506e77315bcSPatrick Sanan           gs = coef*bs*(t0 - ts);
507e77315bcSPatrick Sanan 
508e77315bcSPatrick Sanan           tn = x[k][j+1][i];
509e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
510e77315bcSPatrick Sanan           bn = PetscPowScalar(an,bm1);
511e77315bcSPatrick Sanan           /* dn = bn * an; */
512e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
513e77315bcSPatrick Sanan           gn = coef*bn*(tn - t0);
514e77315bcSPatrick Sanan 
515e77315bcSPatrick Sanan           td = x[k-1][j][i];
516e77315bcSPatrick Sanan           ad = 0.5*(t0 + td);
517e77315bcSPatrick Sanan           bd = PetscPowScalar(ad,bm1);
518e77315bcSPatrick Sanan           /* dd = bd * ad; */
519e77315bcSPatrick Sanan           dd = PetscPowScalar(ad,beta);
520e77315bcSPatrick Sanan           gd = coef*bd*(t0 - td);
521e77315bcSPatrick Sanan 
522e77315bcSPatrick Sanan           tu = x[k+1][j][i];
523e77315bcSPatrick Sanan           au = 0.5*(t0 + tu);
524e77315bcSPatrick Sanan           bu = PetscPowScalar(au,bm1);
525e77315bcSPatrick Sanan           /* du = bu * au; */
526e77315bcSPatrick Sanan           du = PetscPowScalar(au,beta);
527e77315bcSPatrick Sanan           gu = coef*bu*(tu - t0);
528e77315bcSPatrick Sanan 
529e77315bcSPatrick Sanan           c[0].k = k-1; c[0].j = j; c[0].i = i; v[0]   = -hxhydhz*(dd - gd);
530e77315bcSPatrick Sanan           c[1].k = k; c[1].j = j-1; c[1].i = i;
531e77315bcSPatrick Sanan           v[1]   = -hzhxdhy*(ds - gs);
532e77315bcSPatrick Sanan           c[2].k = k; c[2].j = j; c[2].i = i-1;
533e77315bcSPatrick Sanan           v[2]   = -hyhzdhx*(dw - gw);
534e77315bcSPatrick Sanan           c[3].k = k; c[3].j = j; c[3].i = i;
535e77315bcSPatrick Sanan           v[3]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
536e77315bcSPatrick Sanan           c[4].k = k; c[4].j = j; c[4].i = i+1;
537e77315bcSPatrick Sanan           v[4]   = -hyhzdhx*(de + ge);
538e77315bcSPatrick Sanan           c[5].k = k; c[5].j = j+1; c[5].i = i;
539e77315bcSPatrick Sanan           v[5]   = -hzhxdhy*(dn + gn);
540e77315bcSPatrick Sanan           c[6].k = k+1; c[6].j = j; c[6].i = i;
541e77315bcSPatrick Sanan           v[6]   = -hxhydhz*(du + gu);
5425f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(jac,1,&row,7,c,v,INSERT_VALUES));
543e77315bcSPatrick Sanan 
544e77315bcSPatrick Sanan         } else if (i == 0) {
545e77315bcSPatrick Sanan 
546e77315bcSPatrick Sanan           /* left-hand plane boundary */
547e77315bcSPatrick Sanan           tw = tleft;
548e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
549e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
550e77315bcSPatrick Sanan           /* dw = bw * aw */
551e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
552e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
553e77315bcSPatrick Sanan 
554e77315bcSPatrick Sanan           te = x[k][j][i+1];
555e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
556e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
557e77315bcSPatrick Sanan           /* de = be * ae; */
558e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
559e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
560e77315bcSPatrick Sanan 
561e77315bcSPatrick Sanan           /* left-hand bottom edge */
562e77315bcSPatrick Sanan           if (j == 0) {
563e77315bcSPatrick Sanan 
564e77315bcSPatrick Sanan             tn = x[k][j+1][i];
565e77315bcSPatrick Sanan             an = 0.5*(t0 + tn);
566e77315bcSPatrick Sanan             bn = PetscPowScalar(an,bm1);
567e77315bcSPatrick Sanan             /* dn = bn * an; */
568e77315bcSPatrick Sanan             dn = PetscPowScalar(an,beta);
569e77315bcSPatrick Sanan             gn = coef*bn*(tn - t0);
570e77315bcSPatrick Sanan 
571e77315bcSPatrick Sanan             /* left-hand bottom down corner */
572e77315bcSPatrick Sanan             if (k == 0) {
573e77315bcSPatrick Sanan 
574e77315bcSPatrick Sanan               tu = x[k+1][j][i];
575e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
576e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
577e77315bcSPatrick Sanan               /* du = bu * au; */
578e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
579e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
580e77315bcSPatrick Sanan 
581e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j; c[0].i = i;
582e77315bcSPatrick Sanan               v[0]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
583e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i+1;
584e77315bcSPatrick Sanan               v[1]   = -hyhzdhx*(de + ge);
585e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j+1; c[2].i = i;
586e77315bcSPatrick Sanan               v[2]   = -hzhxdhy*(dn + gn);
587e77315bcSPatrick Sanan               c[3].k = k+1; c[3].j = j; c[3].i = i;
588e77315bcSPatrick Sanan               v[3]   = -hxhydhz*(du + gu);
5895f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
590e77315bcSPatrick Sanan 
591e77315bcSPatrick Sanan               /* left-hand bottom interior edge */
592e77315bcSPatrick Sanan             } else if (k < mz-1) {
593e77315bcSPatrick Sanan 
594e77315bcSPatrick Sanan               tu = x[k+1][j][i];
595e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
596e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
597e77315bcSPatrick Sanan               /* du = bu * au; */
598e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
599e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
600e77315bcSPatrick Sanan 
601e77315bcSPatrick Sanan               td = x[k-1][j][i];
602e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
603e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
604e77315bcSPatrick Sanan               /* dd = bd * ad; */
605e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
606e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
607e77315bcSPatrick Sanan 
608e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
609e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
610e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i;
611e77315bcSPatrick Sanan               v[1]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
612e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i+1;
613e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(de + ge);
614e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
615e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
616e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
617e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
6185f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
619e77315bcSPatrick Sanan 
620e77315bcSPatrick Sanan               /* left-hand bottom up corner */
621e77315bcSPatrick Sanan             } else {
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 
630e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
631e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
632e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i;
633e77315bcSPatrick Sanan               v[1]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
634e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i+1;
635e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(de + ge);
636e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
637e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
6385f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
639e77315bcSPatrick Sanan             }
640e77315bcSPatrick Sanan 
641e77315bcSPatrick Sanan             /* left-hand top edge */
642e77315bcSPatrick Sanan           } else if (j == my-1) {
643e77315bcSPatrick Sanan 
644e77315bcSPatrick Sanan             ts = x[k][j-1][i];
645e77315bcSPatrick Sanan             as = 0.5*(t0 + ts);
646e77315bcSPatrick Sanan             bs = PetscPowScalar(as,bm1);
647e77315bcSPatrick Sanan             /* ds = bs * as; */
648e77315bcSPatrick Sanan             ds = PetscPowScalar(as,beta);
649e77315bcSPatrick Sanan             gs = coef*bs*(ts - t0);
650e77315bcSPatrick Sanan 
651e77315bcSPatrick Sanan             /* left-hand top down corner */
652e77315bcSPatrick Sanan             if (k == 0) {
653e77315bcSPatrick Sanan 
654e77315bcSPatrick Sanan               tu = x[k+1][j][i];
655e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
656e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
657e77315bcSPatrick Sanan               /* du = bu * au; */
658e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
659e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
660e77315bcSPatrick Sanan 
661e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j-1; c[0].i = i;
662e77315bcSPatrick Sanan               v[0]   = -hzhxdhy*(ds - gs);
663e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i;
664e77315bcSPatrick Sanan               v[1]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
665e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i+1;
666e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(de + ge);
667e77315bcSPatrick Sanan               c[3].k = k+1; c[3].j = j; c[3].i = i;
668e77315bcSPatrick Sanan               v[3]   = -hxhydhz*(du + gu);
6695f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
670e77315bcSPatrick Sanan 
671e77315bcSPatrick Sanan               /* left-hand top interior edge */
672e77315bcSPatrick Sanan             } else if (k < mz-1) {
673e77315bcSPatrick Sanan 
674e77315bcSPatrick Sanan               tu = x[k+1][j][i];
675e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
676e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
677e77315bcSPatrick Sanan               /* du = bu * au; */
678e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
679e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
680e77315bcSPatrick Sanan 
681e77315bcSPatrick Sanan               td = x[k-1][j][i];
682e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
683e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
684e77315bcSPatrick Sanan               /* dd = bd * ad; */
685e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
686e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
687e77315bcSPatrick Sanan 
688e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
689e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
690e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
691e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
692e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
693e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
694e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i+1;
695e77315bcSPatrick Sanan               v[3]   = -hyhzdhx*(de + ge);
696e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
697e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
6985f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
699e77315bcSPatrick Sanan 
700e77315bcSPatrick Sanan               /* left-hand top up corner */
701e77315bcSPatrick Sanan             } else {
702e77315bcSPatrick Sanan 
703e77315bcSPatrick Sanan               td = x[k-1][j][i];
704e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
705e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
706e77315bcSPatrick Sanan               /* dd = bd * ad; */
707e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
708e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
709e77315bcSPatrick Sanan 
710e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
711e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
712e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
713e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
714e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
715e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
716e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i+1;
717e77315bcSPatrick Sanan               v[3]   = -hyhzdhx*(de + ge);
7185f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
719e77315bcSPatrick Sanan             }
720e77315bcSPatrick Sanan 
721e77315bcSPatrick Sanan           } else {
722e77315bcSPatrick Sanan 
723e77315bcSPatrick Sanan             ts = x[k][j-1][i];
724e77315bcSPatrick Sanan             as = 0.5*(t0 + ts);
725e77315bcSPatrick Sanan             bs = PetscPowScalar(as,bm1);
726e77315bcSPatrick Sanan             /* ds = bs * as; */
727e77315bcSPatrick Sanan             ds = PetscPowScalar(as,beta);
728e77315bcSPatrick Sanan             gs = coef*bs*(t0 - ts);
729e77315bcSPatrick Sanan 
730e77315bcSPatrick Sanan             tn = x[k][j+1][i];
731e77315bcSPatrick Sanan             an = 0.5*(t0 + tn);
732e77315bcSPatrick Sanan             bn = PetscPowScalar(an,bm1);
733e77315bcSPatrick Sanan             /* dn = bn * an; */
734e77315bcSPatrick Sanan             dn = PetscPowScalar(an,beta);
735e77315bcSPatrick Sanan             gn = coef*bn*(tn - t0);
736e77315bcSPatrick Sanan 
737e77315bcSPatrick Sanan             /* left-hand down interior edge */
738e77315bcSPatrick Sanan             if (k == 0) {
739e77315bcSPatrick Sanan 
740e77315bcSPatrick Sanan               tu = x[k+1][j][i];
741e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
742e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
743e77315bcSPatrick Sanan               /* du = bu * au; */
744e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
745e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
746e77315bcSPatrick Sanan 
747e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j-1; c[0].i = i;
748e77315bcSPatrick Sanan               v[0]   = -hzhxdhy*(ds - gs);
749e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i;
750e77315bcSPatrick Sanan               v[1]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
751e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i+1;
752e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(de + ge);
753e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
754e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
755e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
756e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
7575f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
758e77315bcSPatrick Sanan 
759e77315bcSPatrick Sanan             } else if (k == mz-1) { /* left-hand up interior edge */
760e77315bcSPatrick Sanan 
761e77315bcSPatrick Sanan               td = x[k-1][j][i];
762e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
763e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
764e77315bcSPatrick Sanan               /* dd = bd * ad; */
765e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
766e77315bcSPatrick Sanan               gd = coef*bd*(t0 - td);
767e77315bcSPatrick Sanan 
768e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
769e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
770e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
771e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
772e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
773e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
774e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i+1;
775e77315bcSPatrick Sanan               v[3]   = -hyhzdhx*(de + ge);
776e77315bcSPatrick Sanan               c[4].k = k; c[4].j = j+1; c[4].i = i;
777e77315bcSPatrick Sanan               v[4]   = -hzhxdhy*(dn + gn);
7785f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
779e77315bcSPatrick Sanan             } else { /* left-hand interior plane */
780e77315bcSPatrick Sanan 
781e77315bcSPatrick Sanan               td = x[k-1][j][i];
782e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
783e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
784e77315bcSPatrick Sanan               /* dd = bd * ad; */
785e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
786e77315bcSPatrick Sanan               gd = coef*bd*(t0 - td);
787e77315bcSPatrick Sanan 
788e77315bcSPatrick Sanan               tu = x[k+1][j][i];
789e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
790e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
791e77315bcSPatrick Sanan               /* du = bu * au; */
792e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
793e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
794e77315bcSPatrick Sanan 
795e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
796e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
797e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
798e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
799e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
800e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
801e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i+1;
802e77315bcSPatrick Sanan               v[3]   = -hyhzdhx*(de + ge);
803e77315bcSPatrick Sanan               c[4].k = k; c[4].j = j+1; c[4].i = i;
804e77315bcSPatrick Sanan               v[4]   = -hzhxdhy*(dn + gn);
805e77315bcSPatrick Sanan               c[5].k = k+1; c[5].j = j; c[5].i = i;
806e77315bcSPatrick Sanan               v[5]   = -hxhydhz*(du + gu);
8075f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
808e77315bcSPatrick Sanan             }
809e77315bcSPatrick Sanan           }
810e77315bcSPatrick Sanan 
811e77315bcSPatrick Sanan         } else if (i == mx-1) {
812e77315bcSPatrick Sanan 
813e77315bcSPatrick Sanan           /* right-hand plane boundary */
814e77315bcSPatrick Sanan           tw = x[k][j][i-1];
815e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
816e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
817e77315bcSPatrick Sanan           /* dw = bw * aw */
818e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
819e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
820e77315bcSPatrick Sanan 
821e77315bcSPatrick Sanan           te = tright;
822e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
823e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
824e77315bcSPatrick Sanan           /* de = be * ae; */
825e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
826e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
827e77315bcSPatrick Sanan 
828e77315bcSPatrick Sanan           /* right-hand bottom edge */
829e77315bcSPatrick Sanan           if (j == 0) {
830e77315bcSPatrick Sanan 
831e77315bcSPatrick Sanan             tn = x[k][j+1][i];
832e77315bcSPatrick Sanan             an = 0.5*(t0 + tn);
833e77315bcSPatrick Sanan             bn = PetscPowScalar(an,bm1);
834e77315bcSPatrick Sanan             /* dn = bn * an; */
835e77315bcSPatrick Sanan             dn = PetscPowScalar(an,beta);
836e77315bcSPatrick Sanan             gn = coef*bn*(tn - t0);
837e77315bcSPatrick Sanan 
838e77315bcSPatrick Sanan             /* right-hand bottom down corner */
839e77315bcSPatrick Sanan             if (k == 0) {
840e77315bcSPatrick Sanan 
841e77315bcSPatrick Sanan               tu = x[k+1][j][i];
842e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
843e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
844e77315bcSPatrick Sanan               /* du = bu * au; */
845e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
846e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
847e77315bcSPatrick Sanan 
848e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j; c[0].i = i-1;
849e77315bcSPatrick Sanan               v[0]   = -hyhzdhx*(dw - gw);
850e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i;
851e77315bcSPatrick Sanan               v[1]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
852e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j+1; c[2].i = i;
853e77315bcSPatrick Sanan               v[2]   = -hzhxdhy*(dn + gn);
854e77315bcSPatrick Sanan               c[3].k = k+1; c[3].j = j; c[3].i = i;
855e77315bcSPatrick Sanan               v[3]   = -hxhydhz*(du + gu);
8565f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
857e77315bcSPatrick Sanan 
858e77315bcSPatrick Sanan               /* right-hand bottom interior edge */
859e77315bcSPatrick Sanan             } else if (k < mz-1) {
860e77315bcSPatrick Sanan 
861e77315bcSPatrick Sanan               tu = x[k+1][j][i];
862e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
863e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
864e77315bcSPatrick Sanan               /* du = bu * au; */
865e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
866e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
867e77315bcSPatrick Sanan 
868e77315bcSPatrick Sanan               td = x[k-1][j][i];
869e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
870e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
871e77315bcSPatrick Sanan               /* dd = bd * ad; */
872e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
873e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
874e77315bcSPatrick Sanan 
875e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
876e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
877e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i-1;
878e77315bcSPatrick Sanan               v[1]   = -hyhzdhx*(dw - gw);
879e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
880e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
881e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
882e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
883e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
884e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
8855f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
886e77315bcSPatrick Sanan 
887e77315bcSPatrick Sanan               /* right-hand bottom up corner */
888e77315bcSPatrick Sanan             } else {
889e77315bcSPatrick Sanan 
890e77315bcSPatrick Sanan               td = x[k-1][j][i];
891e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
892e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
893e77315bcSPatrick Sanan               /* dd = bd * ad; */
894e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
895e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
896e77315bcSPatrick Sanan 
897e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
898e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
899e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i-1;
900e77315bcSPatrick Sanan               v[1]   = -hyhzdhx*(dw - gw);
901e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
902e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
903e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
904e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
9055f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
906e77315bcSPatrick Sanan             }
907e77315bcSPatrick Sanan 
908e77315bcSPatrick Sanan             /* right-hand top edge */
909e77315bcSPatrick Sanan           } else if (j == my-1) {
910e77315bcSPatrick Sanan 
911e77315bcSPatrick Sanan             ts = x[k][j-1][i];
912e77315bcSPatrick Sanan             as = 0.5*(t0 + ts);
913e77315bcSPatrick Sanan             bs = PetscPowScalar(as,bm1);
914e77315bcSPatrick Sanan             /* ds = bs * as; */
915e77315bcSPatrick Sanan             ds = PetscPowScalar(as,beta);
916e77315bcSPatrick Sanan             gs = coef*bs*(ts - t0);
917e77315bcSPatrick Sanan 
918e77315bcSPatrick Sanan             /* right-hand top down corner */
919e77315bcSPatrick Sanan             if (k == 0) {
920e77315bcSPatrick Sanan 
921e77315bcSPatrick Sanan               tu = x[k+1][j][i];
922e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
923e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
924e77315bcSPatrick Sanan               /* du = bu * au; */
925e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
926e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
927e77315bcSPatrick Sanan 
928e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j-1; c[0].i = i;
929e77315bcSPatrick Sanan               v[0]   = -hzhxdhy*(ds - gs);
930e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i-1;
931e77315bcSPatrick Sanan               v[1]   = -hyhzdhx*(dw - gw);
932e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
933e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
934e77315bcSPatrick Sanan               c[3].k = k+1; c[3].j = j; c[3].i = i;
935e77315bcSPatrick Sanan               v[3]   = -hxhydhz*(du + gu);
9365f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
937e77315bcSPatrick Sanan 
938e77315bcSPatrick Sanan               /* right-hand top interior edge */
939e77315bcSPatrick Sanan             } else if (k < mz-1) {
940e77315bcSPatrick Sanan 
941e77315bcSPatrick Sanan               tu = x[k+1][j][i];
942e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
943e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
944e77315bcSPatrick Sanan               /* du = bu * au; */
945e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
946e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
947e77315bcSPatrick Sanan 
948e77315bcSPatrick Sanan               td = x[k-1][j][i];
949e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
950e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
951e77315bcSPatrick Sanan               /* dd = bd * ad; */
952e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
953e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
954e77315bcSPatrick Sanan 
955e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
956e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
957e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
958e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
959e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i-1;
960e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(dw - gw);
961e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i;
962e77315bcSPatrick Sanan               v[3]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
963e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
964e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
9655f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
966e77315bcSPatrick Sanan 
967e77315bcSPatrick Sanan               /* right-hand top up corner */
968e77315bcSPatrick Sanan             } else {
969e77315bcSPatrick Sanan 
970e77315bcSPatrick Sanan               td = x[k-1][j][i];
971e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
972e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
973e77315bcSPatrick Sanan               /* dd = bd * ad; */
974e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
975e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
976e77315bcSPatrick Sanan 
977e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
978e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
979e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
980e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
981e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i-1;
982e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(dw - gw);
983e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i;
984e77315bcSPatrick Sanan               v[3]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
9855f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
986e77315bcSPatrick Sanan             }
987e77315bcSPatrick Sanan 
988e77315bcSPatrick Sanan           } else {
989e77315bcSPatrick Sanan 
990e77315bcSPatrick Sanan             ts = x[k][j-1][i];
991e77315bcSPatrick Sanan             as = 0.5*(t0 + ts);
992e77315bcSPatrick Sanan             bs = PetscPowScalar(as,bm1);
993e77315bcSPatrick Sanan             /* ds = bs * as; */
994e77315bcSPatrick Sanan             ds = PetscPowScalar(as,beta);
995e77315bcSPatrick Sanan             gs = coef*bs*(t0 - ts);
996e77315bcSPatrick Sanan 
997e77315bcSPatrick Sanan             tn = x[k][j+1][i];
998e77315bcSPatrick Sanan             an = 0.5*(t0 + tn);
999e77315bcSPatrick Sanan             bn = PetscPowScalar(an,bm1);
1000e77315bcSPatrick Sanan             /* dn = bn * an; */
1001e77315bcSPatrick Sanan             dn = PetscPowScalar(an,beta);
1002e77315bcSPatrick Sanan             gn = coef*bn*(tn - t0);
1003e77315bcSPatrick Sanan 
1004e77315bcSPatrick Sanan             /* right-hand down interior edge */
1005e77315bcSPatrick Sanan             if (k == 0) {
1006e77315bcSPatrick Sanan 
1007e77315bcSPatrick Sanan               tu = x[k+1][j][i];
1008e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
1009e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
1010e77315bcSPatrick Sanan               /* du = bu * au; */
1011e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
1012e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
1013e77315bcSPatrick Sanan 
1014e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j-1; c[0].i = i;
1015e77315bcSPatrick Sanan               v[0]   = -hzhxdhy*(ds - gs);
1016e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i-1;
1017e77315bcSPatrick Sanan               v[1]   = -hyhzdhx*(dw - gw);
1018e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
1019e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
1020e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
1021e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
1022e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
1023e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
10245f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1025e77315bcSPatrick Sanan 
1026e77315bcSPatrick Sanan             } else if (k == mz-1) { /* right-hand up interior edge */
1027e77315bcSPatrick Sanan 
1028e77315bcSPatrick Sanan               td = x[k-1][j][i];
1029e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
1030e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
1031e77315bcSPatrick Sanan               /* dd = bd * ad; */
1032e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
1033e77315bcSPatrick Sanan               gd = coef*bd*(t0 - td);
1034e77315bcSPatrick Sanan 
1035e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
1036e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
1037e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
1038e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
1039e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i-1;
1040e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(dw - gw);
1041e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i;
1042e77315bcSPatrick Sanan               v[3]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
1043e77315bcSPatrick Sanan               c[4].k = k; c[4].j = j+1; c[4].i = i;
1044e77315bcSPatrick Sanan               v[4]   = -hzhxdhy*(dn + gn);
10455f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1046e77315bcSPatrick Sanan 
1047e77315bcSPatrick Sanan             } else { /* right-hand interior plane */
1048e77315bcSPatrick Sanan 
1049e77315bcSPatrick Sanan               td = x[k-1][j][i];
1050e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
1051e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
1052e77315bcSPatrick Sanan               /* dd = bd * ad; */
1053e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
1054e77315bcSPatrick Sanan               gd = coef*bd*(t0 - td);
1055e77315bcSPatrick Sanan 
1056e77315bcSPatrick Sanan               tu = x[k+1][j][i];
1057e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
1058e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
1059e77315bcSPatrick Sanan               /* du = bu * au; */
1060e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
1061e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
1062e77315bcSPatrick Sanan 
1063e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
1064e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
1065e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
1066e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
1067e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i-1;
1068e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(dw - gw);
1069e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i;
1070e77315bcSPatrick Sanan               v[3]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
1071e77315bcSPatrick Sanan               c[4].k = k; c[4].j = j+1; c[4].i = i;
1072e77315bcSPatrick Sanan               v[4]   = -hzhxdhy*(dn + gn);
1073e77315bcSPatrick Sanan               c[5].k = k+1; c[5].j = j; c[5].i = i;
1074e77315bcSPatrick Sanan               v[5]   = -hxhydhz*(du + gu);
10755f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
1076e77315bcSPatrick Sanan             }
1077e77315bcSPatrick Sanan           }
1078e77315bcSPatrick Sanan 
1079e77315bcSPatrick Sanan         } else if (j == 0) {
1080e77315bcSPatrick Sanan 
1081e77315bcSPatrick Sanan           tw = x[k][j][i-1];
1082e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
1083e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
1084e77315bcSPatrick Sanan           /* dw = bw * aw */
1085e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
1086e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
1087e77315bcSPatrick Sanan 
1088e77315bcSPatrick Sanan           te = x[k][j][i+1];
1089e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
1090e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
1091e77315bcSPatrick Sanan           /* de = be * ae; */
1092e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
1093e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
1094e77315bcSPatrick Sanan 
1095e77315bcSPatrick Sanan           tn = x[k][j+1][i];
1096e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
1097e77315bcSPatrick Sanan           bn = PetscPowScalar(an,bm1);
1098e77315bcSPatrick Sanan           /* dn = bn * an; */
1099e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
1100e77315bcSPatrick Sanan           gn = coef*bn*(tn - t0);
1101e77315bcSPatrick Sanan 
1102e77315bcSPatrick Sanan           /* bottom down interior edge */
1103e77315bcSPatrick Sanan           if (k == 0) {
1104e77315bcSPatrick Sanan 
1105e77315bcSPatrick Sanan             tu = x[k+1][j][i];
1106e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
1107e77315bcSPatrick Sanan             bu = PetscPowScalar(au,bm1);
1108e77315bcSPatrick Sanan             /* du = bu * au; */
1109e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
1110e77315bcSPatrick Sanan             gu = coef*bu*(tu - t0);
1111e77315bcSPatrick Sanan 
1112e77315bcSPatrick Sanan             c[0].k = k; c[0].j = j; c[0].i = i-1;
1113e77315bcSPatrick Sanan             v[0]   = -hyhzdhx*(dw - gw);
1114e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j; c[1].i = i;
1115e77315bcSPatrick Sanan             v[1]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
1116e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i+1;
1117e77315bcSPatrick Sanan             v[2]   = -hyhzdhx*(de + ge);
1118e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j+1; c[3].i = i;
1119e77315bcSPatrick Sanan             v[3]   = -hzhxdhy*(dn + gn);
1120e77315bcSPatrick Sanan             c[4].k = k+1; c[4].j = j; c[4].i = i;
1121e77315bcSPatrick Sanan             v[4]   = -hxhydhz*(du + gu);
11225f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1123e77315bcSPatrick Sanan 
1124e77315bcSPatrick Sanan           } else if (k == mz-1) { /* bottom up interior edge */
1125e77315bcSPatrick Sanan 
1126e77315bcSPatrick Sanan             td = x[k-1][j][i];
1127e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
1128e77315bcSPatrick Sanan             bd = PetscPowScalar(ad,bm1);
1129e77315bcSPatrick Sanan             /* dd = bd * ad; */
1130e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
1131e77315bcSPatrick Sanan             gd = coef*bd*(td - t0);
1132e77315bcSPatrick Sanan 
1133e77315bcSPatrick Sanan             c[0].k = k-1; c[0].j = j; c[0].i = i;
1134e77315bcSPatrick Sanan             v[0]   = -hxhydhz*(dd - gd);
1135e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j; c[1].i = i-1;
1136e77315bcSPatrick Sanan             v[1]   = -hyhzdhx*(dw - gw);
1137e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i;
1138e77315bcSPatrick Sanan             v[2]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
1139e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j; c[3].i = i+1;
1140e77315bcSPatrick Sanan             v[3]   = -hyhzdhx*(de + ge);
1141e77315bcSPatrick Sanan             c[4].k = k; c[4].j = j+1; c[4].i = i;
1142e77315bcSPatrick Sanan             v[4]   = -hzhxdhy*(dn + gn);
11435f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1144e77315bcSPatrick Sanan 
1145e77315bcSPatrick Sanan           } else { /* bottom interior plane */
1146e77315bcSPatrick Sanan 
1147e77315bcSPatrick Sanan             tu = x[k+1][j][i];
1148e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
1149e77315bcSPatrick Sanan             bu = PetscPowScalar(au,bm1);
1150e77315bcSPatrick Sanan             /* du = bu * au; */
1151e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
1152e77315bcSPatrick Sanan             gu = coef*bu*(tu - t0);
1153e77315bcSPatrick Sanan 
1154e77315bcSPatrick Sanan             td = x[k-1][j][i];
1155e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
1156e77315bcSPatrick Sanan             bd = PetscPowScalar(ad,bm1);
1157e77315bcSPatrick Sanan             /* dd = bd * ad; */
1158e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
1159e77315bcSPatrick Sanan             gd = coef*bd*(td - t0);
1160e77315bcSPatrick Sanan 
1161e77315bcSPatrick Sanan             c[0].k = k-1; c[0].j = j; c[0].i = i;
1162e77315bcSPatrick Sanan             v[0]   = -hxhydhz*(dd - gd);
1163e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j; c[1].i = i-1;
1164e77315bcSPatrick Sanan             v[1]   = -hyhzdhx*(dw - gw);
1165e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i;
1166e77315bcSPatrick Sanan             v[2]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
1167e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j; c[3].i = i+1;
1168e77315bcSPatrick Sanan             v[3]   = -hyhzdhx*(de + ge);
1169e77315bcSPatrick Sanan             c[4].k = k; c[4].j = j+1; c[4].i = i;
1170e77315bcSPatrick Sanan             v[4]   = -hzhxdhy*(dn + gn);
1171e77315bcSPatrick Sanan             c[5].k = k+1; c[5].j = j; c[5].i = i;
1172e77315bcSPatrick Sanan             v[5]   = -hxhydhz*(du + gu);
11735f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
1174e77315bcSPatrick Sanan           }
1175e77315bcSPatrick Sanan 
1176e77315bcSPatrick Sanan         } else if (j == my-1) {
1177e77315bcSPatrick Sanan 
1178e77315bcSPatrick Sanan           tw = x[k][j][i-1];
1179e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
1180e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
1181e77315bcSPatrick Sanan           /* dw = bw * aw */
1182e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
1183e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
1184e77315bcSPatrick Sanan 
1185e77315bcSPatrick Sanan           te = x[k][j][i+1];
1186e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
1187e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
1188e77315bcSPatrick Sanan           /* de = be * ae; */
1189e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
1190e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
1191e77315bcSPatrick Sanan 
1192e77315bcSPatrick Sanan           ts = x[k][j-1][i];
1193e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
1194e77315bcSPatrick Sanan           bs = PetscPowScalar(as,bm1);
1195e77315bcSPatrick Sanan           /* ds = bs * as; */
1196e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
1197e77315bcSPatrick Sanan           gs = coef*bs*(t0 - ts);
1198e77315bcSPatrick Sanan 
1199e77315bcSPatrick Sanan           /* top down interior edge */
1200e77315bcSPatrick Sanan           if (k == 0) {
1201e77315bcSPatrick Sanan 
1202e77315bcSPatrick Sanan             tu = x[k+1][j][i];
1203e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
1204e77315bcSPatrick Sanan             bu = PetscPowScalar(au,bm1);
1205e77315bcSPatrick Sanan             /* du = bu * au; */
1206e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
1207e77315bcSPatrick Sanan             gu = coef*bu*(tu - t0);
1208e77315bcSPatrick Sanan 
1209e77315bcSPatrick Sanan             c[0].k = k; c[0].j = j-1; c[0].i = i;
1210e77315bcSPatrick Sanan             v[0]   = -hzhxdhy*(ds - gs);
1211e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j; c[1].i = i-1;
1212e77315bcSPatrick Sanan             v[1]   = -hyhzdhx*(dw - gw);
1213e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i;
1214e77315bcSPatrick Sanan             v[2]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
1215e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j; c[3].i = i+1;
1216e77315bcSPatrick Sanan             v[3]   = -hyhzdhx*(de + ge);
1217e77315bcSPatrick Sanan             c[4].k = k+1; c[4].j = j; c[4].i = i;
1218e77315bcSPatrick Sanan             v[4]   = -hxhydhz*(du + gu);
12195f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1220e77315bcSPatrick Sanan 
1221e77315bcSPatrick Sanan           } else if (k == mz-1) { /* top up interior edge */
1222e77315bcSPatrick Sanan 
1223e77315bcSPatrick Sanan             td = x[k-1][j][i];
1224e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
1225e77315bcSPatrick Sanan             bd = PetscPowScalar(ad,bm1);
1226e77315bcSPatrick Sanan             /* dd = bd * ad; */
1227e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
1228e77315bcSPatrick Sanan             gd = coef*bd*(td - t0);
1229e77315bcSPatrick Sanan 
1230e77315bcSPatrick Sanan             c[0].k = k-1; c[0].j = j; c[0].i = i;
1231e77315bcSPatrick Sanan             v[0]   = -hxhydhz*(dd - gd);
1232e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j-1; c[1].i = i;
1233e77315bcSPatrick Sanan             v[1]   = -hzhxdhy*(ds - gs);
1234e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i-1;
1235e77315bcSPatrick Sanan             v[2]   = -hyhzdhx*(dw - gw);
1236e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j; c[3].i = i;
1237e77315bcSPatrick Sanan             v[3]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
1238e77315bcSPatrick Sanan             c[4].k = k; c[4].j = j; c[4].i = i+1;
1239e77315bcSPatrick Sanan             v[4]   = -hyhzdhx*(de + ge);
12405f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1241e77315bcSPatrick Sanan 
1242e77315bcSPatrick Sanan           } else { /* top interior plane */
1243e77315bcSPatrick Sanan 
1244e77315bcSPatrick Sanan             tu = x[k+1][j][i];
1245e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
1246e77315bcSPatrick Sanan             bu = PetscPowScalar(au,bm1);
1247e77315bcSPatrick Sanan             /* du = bu * au; */
1248e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
1249e77315bcSPatrick Sanan             gu = coef*bu*(tu - t0);
1250e77315bcSPatrick Sanan 
1251e77315bcSPatrick Sanan             td = x[k-1][j][i];
1252e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
1253e77315bcSPatrick Sanan             bd = PetscPowScalar(ad,bm1);
1254e77315bcSPatrick Sanan             /* dd = bd * ad; */
1255e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
1256e77315bcSPatrick Sanan             gd = coef*bd*(td - t0);
1257e77315bcSPatrick Sanan 
1258e77315bcSPatrick Sanan             c[0].k = k-1; c[0].j = j; c[0].i = i;
1259e77315bcSPatrick Sanan             v[0]   = -hxhydhz*(dd - gd);
1260e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j-1; c[1].i = i;
1261e77315bcSPatrick Sanan             v[1]   = -hzhxdhy*(ds - gs);
1262e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i-1;
1263e77315bcSPatrick Sanan             v[2]   = -hyhzdhx*(dw - gw);
1264e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j; c[3].i = i;
1265e77315bcSPatrick Sanan             v[3]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
1266e77315bcSPatrick Sanan             c[4].k = k; c[4].j = j; c[4].i = i+1;
1267e77315bcSPatrick Sanan             v[4]   = -hyhzdhx*(de + ge);
1268e77315bcSPatrick Sanan             c[5].k = k+1; c[5].j = j; c[5].i = i;
1269e77315bcSPatrick Sanan             v[5]   = -hxhydhz*(du + gu);
12705f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
1271e77315bcSPatrick Sanan           }
1272e77315bcSPatrick Sanan 
1273e77315bcSPatrick Sanan         } else if (k == 0) {
1274e77315bcSPatrick Sanan 
1275e77315bcSPatrick Sanan           /* down interior plane */
1276e77315bcSPatrick Sanan 
1277e77315bcSPatrick Sanan           tw = x[k][j][i-1];
1278e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
1279e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
1280e77315bcSPatrick Sanan           /* dw = bw * aw */
1281e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
1282e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
1283e77315bcSPatrick Sanan 
1284e77315bcSPatrick Sanan           te = x[k][j][i+1];
1285e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
1286e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
1287e77315bcSPatrick Sanan           /* de = be * ae; */
1288e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
1289e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
1290e77315bcSPatrick Sanan 
1291e77315bcSPatrick Sanan           ts = x[k][j-1][i];
1292e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
1293e77315bcSPatrick Sanan           bs = PetscPowScalar(as,bm1);
1294e77315bcSPatrick Sanan           /* ds = bs * as; */
1295e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
1296e77315bcSPatrick Sanan           gs = coef*bs*(t0 - ts);
1297e77315bcSPatrick Sanan 
1298e77315bcSPatrick Sanan           tn = x[k][j+1][i];
1299e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
1300e77315bcSPatrick Sanan           bn = PetscPowScalar(an,bm1);
1301e77315bcSPatrick Sanan           /* dn = bn * an; */
1302e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
1303e77315bcSPatrick Sanan           gn = coef*bn*(tn - t0);
1304e77315bcSPatrick Sanan 
1305e77315bcSPatrick Sanan           tu = x[k+1][j][i];
1306e77315bcSPatrick Sanan           au = 0.5*(t0 + tu);
1307e77315bcSPatrick Sanan           bu = PetscPowScalar(au,bm1);
1308e77315bcSPatrick Sanan           /* du = bu * au; */
1309e77315bcSPatrick Sanan           du = PetscPowScalar(au,beta);
1310e77315bcSPatrick Sanan           gu = coef*bu*(tu - t0);
1311e77315bcSPatrick Sanan 
1312e77315bcSPatrick Sanan           c[0].k = k; c[0].j = j-1; c[0].i = i;
1313e77315bcSPatrick Sanan           v[0]   = -hzhxdhy*(ds - gs);
1314e77315bcSPatrick Sanan           c[1].k = k; c[1].j = j; c[1].i = i-1;
1315e77315bcSPatrick Sanan           v[1]   = -hyhzdhx*(dw - gw);
1316e77315bcSPatrick Sanan           c[2].k = k; c[2].j = j; c[2].i = i;
1317e77315bcSPatrick Sanan           v[2]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
1318e77315bcSPatrick Sanan           c[3].k = k; c[3].j = j; c[3].i = i+1;
1319e77315bcSPatrick Sanan           v[3]   = -hyhzdhx*(de + ge);
1320e77315bcSPatrick Sanan           c[4].k = k; c[4].j = j+1; c[4].i = i;
1321e77315bcSPatrick Sanan           v[4]   = -hzhxdhy*(dn + gn);
1322e77315bcSPatrick Sanan           c[5].k = k+1; c[5].j = j; c[5].i = i;
1323e77315bcSPatrick Sanan           v[5]   = -hxhydhz*(du + gu);
13245f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
1325e77315bcSPatrick Sanan 
1326e77315bcSPatrick Sanan         } else if (k == mz-1) {
1327e77315bcSPatrick Sanan 
1328e77315bcSPatrick Sanan           /* up interior plane */
1329e77315bcSPatrick Sanan 
1330e77315bcSPatrick Sanan           tw = x[k][j][i-1];
1331e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
1332e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
1333e77315bcSPatrick Sanan           /* dw = bw * aw */
1334e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
1335e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
1336e77315bcSPatrick Sanan 
1337e77315bcSPatrick Sanan           te = x[k][j][i+1];
1338e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
1339e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
1340e77315bcSPatrick Sanan           /* de = be * ae; */
1341e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
1342e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
1343e77315bcSPatrick Sanan 
1344e77315bcSPatrick Sanan           ts = x[k][j-1][i];
1345e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
1346e77315bcSPatrick Sanan           bs = PetscPowScalar(as,bm1);
1347e77315bcSPatrick Sanan           /* ds = bs * as; */
1348e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
1349e77315bcSPatrick Sanan           gs = coef*bs*(t0 - ts);
1350e77315bcSPatrick Sanan 
1351e77315bcSPatrick Sanan           tn = x[k][j+1][i];
1352e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
1353e77315bcSPatrick Sanan           bn = PetscPowScalar(an,bm1);
1354e77315bcSPatrick Sanan           /* dn = bn * an; */
1355e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
1356e77315bcSPatrick Sanan           gn = coef*bn*(tn - t0);
1357e77315bcSPatrick Sanan 
1358e77315bcSPatrick Sanan           td = x[k-1][j][i];
1359e77315bcSPatrick Sanan           ad = 0.5*(t0 + td);
1360e77315bcSPatrick Sanan           bd = PetscPowScalar(ad,bm1);
1361e77315bcSPatrick Sanan           /* dd = bd * ad; */
1362e77315bcSPatrick Sanan           dd = PetscPowScalar(ad,beta);
1363e77315bcSPatrick Sanan           gd = coef*bd*(t0 - td);
1364e77315bcSPatrick Sanan 
1365e77315bcSPatrick Sanan           c[0].k = k-1; c[0].j = j; c[0].i = i;
1366e77315bcSPatrick Sanan           v[0]   = -hxhydhz*(dd - gd);
1367e77315bcSPatrick Sanan           c[1].k = k; c[1].j = j-1; c[1].i = i;
1368e77315bcSPatrick Sanan           v[1]   = -hzhxdhy*(ds - gs);
1369e77315bcSPatrick Sanan           c[2].k = k; c[2].j = j; c[2].i = i-1;
1370e77315bcSPatrick Sanan           v[2]   = -hyhzdhx*(dw - gw);
1371e77315bcSPatrick Sanan           c[3].k = k; c[3].j = j; c[3].i = i;
1372e77315bcSPatrick Sanan           v[3]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
1373e77315bcSPatrick Sanan           c[4].k = k; c[4].j = j; c[4].i = i+1;
1374e77315bcSPatrick Sanan           v[4]   = -hyhzdhx*(de + ge);
1375e77315bcSPatrick Sanan           c[5].k = k; c[5].j = j+1; c[5].i = i;
1376e77315bcSPatrick Sanan           v[5]   = -hzhxdhy*(dn + gn);
13775f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
1378e77315bcSPatrick Sanan         }
1379e77315bcSPatrick Sanan       }
1380e77315bcSPatrick Sanan     }
1381e77315bcSPatrick Sanan   }
13825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
13835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
1384e77315bcSPatrick Sanan   if (jac != J) {
13855f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
13865f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
1387e77315bcSPatrick Sanan   }
13885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,localX,&x));
13895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localX));
1390e77315bcSPatrick Sanan 
13915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops((41.0 + 8.0*POWFLOP)*xm*ym));
1392e77315bcSPatrick Sanan   PetscFunctionReturn(0);
1393e77315bcSPatrick Sanan }
1394e77315bcSPatrick Sanan 
1395e77315bcSPatrick Sanan /*TEST
1396e77315bcSPatrick Sanan 
1397e77315bcSPatrick Sanan    test:
1398e77315bcSPatrick Sanan       nsize: 4
1399e77315bcSPatrick 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
1400e77315bcSPatrick Sanan       requires: !single
1401e77315bcSPatrick Sanan 
1402e77315bcSPatrick Sanan TEST*/
1403