xref: /petsc/src/snes/tests/ex20.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
64e77315bcSPatrick Sanan   PetscInt       its,lits;
65e77315bcSPatrick Sanan   PetscReal      litspit;
66e77315bcSPatrick Sanan   DM             da;
67e77315bcSPatrick Sanan 
68e77315bcSPatrick Sanan   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
69e77315bcSPatrick Sanan   /* set problem parameters */
70e77315bcSPatrick Sanan   user.tleft  = 1.0;
71e77315bcSPatrick Sanan   user.tright = 0.1;
72e77315bcSPatrick Sanan   user.beta   = 2.5;
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tleft",&user.tleft,NULL));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tright",&user.tright,NULL));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-beta",&user.beta,NULL));
76e77315bcSPatrick Sanan   user.bm1    = user.beta - 1.0;
77e77315bcSPatrick Sanan   user.coef   = user.beta/2.0;
78e77315bcSPatrick Sanan 
79e77315bcSPatrick Sanan   /*
80e77315bcSPatrick Sanan       Set the DMDA (grid structure) for the grids.
81e77315bcSPatrick Sanan   */
82*5f80ce2aSJacob 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));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(da,&user));
86e77315bcSPatrick Sanan 
87e77315bcSPatrick Sanan   /*
88e77315bcSPatrick Sanan      Create the nonlinear solver
89e77315bcSPatrick Sanan   */
90*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes,da));
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,NULL,FormFunction,&user));
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,NULL,NULL,FormJacobian,&user));
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
95*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetComputeInitialGuess(snes,FormInitialGuess,NULL));
96e77315bcSPatrick Sanan 
97*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,NULL));
98*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetLinearSolveIterations(snes,&lits));
100e77315bcSPatrick Sanan   litspit = ((PetscReal)lits)/((PetscReal)its);
101*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its));
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of Linear iterations = %D\n",lits));
103*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Average Linear its / SNES = %e\n",(double)litspit));
104e77315bcSPatrick Sanan 
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
107e77315bcSPatrick Sanan   ierr = PetscFinalize();
108e77315bcSPatrick Sanan   return ierr;
109e77315bcSPatrick Sanan }
110e77315bcSPatrick Sanan /* --------------------  Form initial approximation ----------------- */
111e77315bcSPatrick Sanan PetscErrorCode FormInitialGuess(SNES snes,Vec X,void *ctx)
112e77315bcSPatrick Sanan {
113e77315bcSPatrick Sanan   AppCtx         *user;
114e77315bcSPatrick Sanan   PetscInt       i,j,k,xs,ys,xm,ym,zs,zm;
115e77315bcSPatrick Sanan   PetscScalar    ***x;
116e77315bcSPatrick Sanan   DM             da;
117e77315bcSPatrick Sanan 
118e77315bcSPatrick Sanan   PetscFunctionBeginUser;
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(da,&user));
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
123e77315bcSPatrick Sanan 
124e77315bcSPatrick Sanan   /* Compute initial guess */
125e77315bcSPatrick Sanan   for (k=zs; k<zs+zm; k++) {
126e77315bcSPatrick Sanan     for (j=ys; j<ys+ym; j++) {
127e77315bcSPatrick Sanan       for (i=xs; i<xs+xm; i++) {
128e77315bcSPatrick Sanan         x[k][j][i] = user->tleft;
129e77315bcSPatrick Sanan       }
130e77315bcSPatrick Sanan     }
131e77315bcSPatrick Sanan   }
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
133e77315bcSPatrick Sanan   PetscFunctionReturn(0);
134e77315bcSPatrick Sanan }
135e77315bcSPatrick Sanan /* --------------------  Evaluate Function F(x) --------------------- */
136e77315bcSPatrick Sanan PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
137e77315bcSPatrick Sanan {
138e77315bcSPatrick Sanan   AppCtx         *user = (AppCtx*)ptr;
139e77315bcSPatrick Sanan   PetscInt       i,j,k,mx,my,mz,xs,ys,zs,xm,ym,zm;
140e77315bcSPatrick Sanan   PetscScalar    zero = 0.0,one = 1.0;
141e77315bcSPatrick Sanan   PetscScalar    hx,hy,hz,hxhydhz,hyhzdhx,hzhxdhy;
142e77315bcSPatrick 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;
143e77315bcSPatrick Sanan   PetscScalar    tleft,tright,beta,td,ad,dd,fd=0.0,tu,au,du=0.0,fu=0.0;
144e77315bcSPatrick Sanan   PetscScalar    ***x,***f;
145e77315bcSPatrick Sanan   Vec            localX;
146e77315bcSPatrick Sanan   DM             da;
147e77315bcSPatrick Sanan 
148e77315bcSPatrick Sanan   PetscFunctionBeginUser;
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
150*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localX));
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
152e77315bcSPatrick Sanan   hx      = one/(PetscReal)(mx-1);  hy    = one/(PetscReal)(my-1);  hz = one/(PetscReal)(mz-1);
153e77315bcSPatrick Sanan   hxhydhz = hx*hy/hz;   hyhzdhx = hy*hz/hx;   hzhxdhy = hz*hx/hy;
154e77315bcSPatrick Sanan   tleft   = user->tleft;         tright = user->tright;
155e77315bcSPatrick Sanan   beta    = user->beta;
156e77315bcSPatrick Sanan 
157e77315bcSPatrick Sanan   /* Get ghost points */
158*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,localX,&x));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,F,&f));
163e77315bcSPatrick Sanan 
164e77315bcSPatrick Sanan   /* Evaluate function */
165e77315bcSPatrick Sanan   for (k=zs; k<zs+zm; k++) {
166e77315bcSPatrick Sanan     for (j=ys; j<ys+ym; j++) {
167e77315bcSPatrick Sanan       for (i=xs; i<xs+xm; i++) {
168e77315bcSPatrick Sanan         t0 = x[k][j][i];
169e77315bcSPatrick Sanan 
170e77315bcSPatrick Sanan         if (i > 0 && i < mx-1 && j > 0 && j < my-1 && k > 0 && k < mz-1) {
171e77315bcSPatrick Sanan 
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 
206e77315bcSPatrick Sanan           /* left-hand (west) boundary */
207e77315bcSPatrick Sanan           tw = tleft;
208e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
209e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
210e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
211e77315bcSPatrick Sanan 
212e77315bcSPatrick Sanan           te = x[k][j][i+1];
213e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
214e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
215e77315bcSPatrick Sanan           fe = de*(te - t0);
216e77315bcSPatrick Sanan 
217e77315bcSPatrick Sanan           if (j > 0) {
218e77315bcSPatrick Sanan             ts = x[k][j-1][i];
219e77315bcSPatrick Sanan             as = 0.5*(t0 + ts);
220e77315bcSPatrick Sanan             ds = PetscPowScalar(as,beta);
221e77315bcSPatrick Sanan             fs = ds*(t0 - ts);
222e77315bcSPatrick Sanan           } else {
223e77315bcSPatrick Sanan             fs = zero;
224e77315bcSPatrick Sanan           }
225e77315bcSPatrick Sanan 
226e77315bcSPatrick Sanan           if (j < my-1) {
227e77315bcSPatrick Sanan             tn = x[k][j+1][i];
228e77315bcSPatrick Sanan             an = 0.5*(t0 + tn);
229e77315bcSPatrick Sanan             dn = PetscPowScalar(an,beta);
230e77315bcSPatrick Sanan             fn = dn*(tn - t0);
231e77315bcSPatrick Sanan           } else {
232e77315bcSPatrick Sanan             fn = zero;
233e77315bcSPatrick Sanan           }
234e77315bcSPatrick Sanan 
235e77315bcSPatrick Sanan           if (k > 0) {
236e77315bcSPatrick Sanan             td = x[k-1][j][i];
237e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
238e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
239e77315bcSPatrick Sanan             fd = dd*(t0 - td);
240e77315bcSPatrick Sanan           } else {
241e77315bcSPatrick Sanan             fd = zero;
242e77315bcSPatrick Sanan           }
243e77315bcSPatrick Sanan 
244e77315bcSPatrick Sanan           if (k < mz-1) {
245e77315bcSPatrick Sanan             tu = x[k+1][j][i];
246e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
247e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
248e77315bcSPatrick Sanan             fu = du*(tu - t0);
249e77315bcSPatrick Sanan           } else {
250e77315bcSPatrick Sanan             fu = zero;
251e77315bcSPatrick Sanan           }
252e77315bcSPatrick Sanan 
253e77315bcSPatrick Sanan         } else if (i == mx-1) {
254e77315bcSPatrick Sanan 
255e77315bcSPatrick Sanan           /* right-hand (east) boundary */
256e77315bcSPatrick Sanan           tw = x[k][j][i-1];
257e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
258e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
259e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
260e77315bcSPatrick Sanan 
261e77315bcSPatrick Sanan           te = tright;
262e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
263e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
264e77315bcSPatrick Sanan           fe = de*(te - t0);
265e77315bcSPatrick Sanan 
266e77315bcSPatrick Sanan           if (j > 0) {
267e77315bcSPatrick Sanan             ts = x[k][j-1][i];
268e77315bcSPatrick Sanan             as = 0.5*(t0 + ts);
269e77315bcSPatrick Sanan             ds = PetscPowScalar(as,beta);
270e77315bcSPatrick Sanan             fs = ds*(t0 - ts);
271e77315bcSPatrick Sanan           } else {
272e77315bcSPatrick Sanan             fs = zero;
273e77315bcSPatrick Sanan           }
274e77315bcSPatrick Sanan 
275e77315bcSPatrick Sanan           if (j < my-1) {
276e77315bcSPatrick Sanan             tn = x[k][j+1][i];
277e77315bcSPatrick Sanan             an = 0.5*(t0 + tn);
278e77315bcSPatrick Sanan             dn = PetscPowScalar(an,beta);
279e77315bcSPatrick Sanan             fn = dn*(tn - t0);
280e77315bcSPatrick Sanan           } else {
281e77315bcSPatrick Sanan             fn = zero;
282e77315bcSPatrick Sanan           }
283e77315bcSPatrick Sanan 
284e77315bcSPatrick Sanan           if (k > 0) {
285e77315bcSPatrick Sanan             td = x[k-1][j][i];
286e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
287e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
288e77315bcSPatrick Sanan             fd = dd*(t0 - td);
289e77315bcSPatrick Sanan           } else {
290e77315bcSPatrick Sanan             fd = zero;
291e77315bcSPatrick Sanan           }
292e77315bcSPatrick Sanan 
293e77315bcSPatrick Sanan           if (k < mz-1) {
294e77315bcSPatrick Sanan             tu = x[k+1][j][i];
295e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
296e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
297e77315bcSPatrick Sanan             fu = du*(tu - t0);
298e77315bcSPatrick Sanan           } else {
299e77315bcSPatrick Sanan             fu = zero;
300e77315bcSPatrick Sanan           }
301e77315bcSPatrick Sanan 
302e77315bcSPatrick Sanan         } else if (j == 0) {
303e77315bcSPatrick Sanan 
304e77315bcSPatrick Sanan           /* bottom (south) boundary, and i <> 0 or mx-1 */
305e77315bcSPatrick Sanan           tw = x[k][j][i-1];
306e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
307e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
308e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
309e77315bcSPatrick Sanan 
310e77315bcSPatrick Sanan           te = x[k][j][i+1];
311e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
312e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
313e77315bcSPatrick Sanan           fe = de*(te - t0);
314e77315bcSPatrick Sanan 
315e77315bcSPatrick Sanan           fs = zero;
316e77315bcSPatrick Sanan 
317e77315bcSPatrick Sanan           tn = x[k][j+1][i];
318e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
319e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
320e77315bcSPatrick Sanan           fn = dn*(tn - t0);
321e77315bcSPatrick Sanan 
322e77315bcSPatrick Sanan           if (k > 0) {
323e77315bcSPatrick Sanan             td = x[k-1][j][i];
324e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
325e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
326e77315bcSPatrick Sanan             fd = dd*(t0 - td);
327e77315bcSPatrick Sanan           } else {
328e77315bcSPatrick Sanan             fd = zero;
329e77315bcSPatrick Sanan           }
330e77315bcSPatrick Sanan 
331e77315bcSPatrick Sanan           if (k < mz-1) {
332e77315bcSPatrick Sanan             tu = x[k+1][j][i];
333e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
334e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
335e77315bcSPatrick Sanan             fu = du*(tu - t0);
336e77315bcSPatrick Sanan           } else {
337e77315bcSPatrick Sanan             fu = zero;
338e77315bcSPatrick Sanan           }
339e77315bcSPatrick Sanan 
340e77315bcSPatrick Sanan         } else if (j == my-1) {
341e77315bcSPatrick Sanan 
342e77315bcSPatrick Sanan           /* top (north) boundary, and i <> 0 or mx-1 */
343e77315bcSPatrick Sanan           tw = x[k][j][i-1];
344e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
345e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
346e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
347e77315bcSPatrick Sanan 
348e77315bcSPatrick Sanan           te = x[k][j][i+1];
349e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
350e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
351e77315bcSPatrick Sanan           fe = de*(te - t0);
352e77315bcSPatrick Sanan 
353e77315bcSPatrick Sanan           ts = x[k][j-1][i];
354e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
355e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
356e77315bcSPatrick Sanan           fs = ds*(t0 - ts);
357e77315bcSPatrick Sanan 
358e77315bcSPatrick Sanan           fn = zero;
359e77315bcSPatrick Sanan 
360e77315bcSPatrick Sanan           if (k > 0) {
361e77315bcSPatrick Sanan             td = x[k-1][j][i];
362e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
363e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
364e77315bcSPatrick Sanan             fd = dd*(t0 - td);
365e77315bcSPatrick Sanan           } else {
366e77315bcSPatrick Sanan             fd = zero;
367e77315bcSPatrick Sanan           }
368e77315bcSPatrick Sanan 
369e77315bcSPatrick Sanan           if (k < mz-1) {
370e77315bcSPatrick Sanan             tu = x[k+1][j][i];
371e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
372e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
373e77315bcSPatrick Sanan             fu = du*(tu - t0);
374e77315bcSPatrick Sanan           } else {
375e77315bcSPatrick Sanan             fu = zero;
376e77315bcSPatrick Sanan           }
377e77315bcSPatrick Sanan 
378e77315bcSPatrick Sanan         } else if (k == 0) {
379e77315bcSPatrick Sanan 
380e77315bcSPatrick Sanan           /* down boundary (interior only) */
381e77315bcSPatrick Sanan           tw = x[k][j][i-1];
382e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
383e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
384e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
385e77315bcSPatrick Sanan 
386e77315bcSPatrick Sanan           te = x[k][j][i+1];
387e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
388e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
389e77315bcSPatrick Sanan           fe = de*(te - t0);
390e77315bcSPatrick Sanan 
391e77315bcSPatrick Sanan           ts = x[k][j-1][i];
392e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
393e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
394e77315bcSPatrick Sanan           fs = ds*(t0 - ts);
395e77315bcSPatrick Sanan 
396e77315bcSPatrick Sanan           tn = x[k][j+1][i];
397e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
398e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
399e77315bcSPatrick Sanan           fn = dn*(tn - t0);
400e77315bcSPatrick Sanan 
401e77315bcSPatrick Sanan           fd = zero;
402e77315bcSPatrick Sanan 
403e77315bcSPatrick Sanan           tu = x[k+1][j][i];
404e77315bcSPatrick Sanan           au = 0.5*(t0 + tu);
405e77315bcSPatrick Sanan           du = PetscPowScalar(au,beta);
406e77315bcSPatrick Sanan           fu = du*(tu - t0);
407e77315bcSPatrick Sanan 
408e77315bcSPatrick Sanan         } else if (k == mz-1) {
409e77315bcSPatrick Sanan 
410e77315bcSPatrick Sanan           /* up boundary (interior only) */
411e77315bcSPatrick Sanan           tw = x[k][j][i-1];
412e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
413e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
414e77315bcSPatrick Sanan           fw = dw*(t0 - tw);
415e77315bcSPatrick Sanan 
416e77315bcSPatrick Sanan           te = x[k][j][i+1];
417e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
418e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
419e77315bcSPatrick Sanan           fe = de*(te - t0);
420e77315bcSPatrick Sanan 
421e77315bcSPatrick Sanan           ts = x[k][j-1][i];
422e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
423e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
424e77315bcSPatrick Sanan           fs = ds*(t0 - ts);
425e77315bcSPatrick Sanan 
426e77315bcSPatrick Sanan           tn = x[k][j+1][i];
427e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
428e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
429e77315bcSPatrick Sanan           fn = dn*(tn - t0);
430e77315bcSPatrick Sanan 
431e77315bcSPatrick Sanan           td = x[k-1][j][i];
432e77315bcSPatrick Sanan           ad = 0.5*(t0 + td);
433e77315bcSPatrick Sanan           dd = PetscPowScalar(ad,beta);
434e77315bcSPatrick Sanan           fd = dd*(t0 - td);
435e77315bcSPatrick Sanan 
436e77315bcSPatrick Sanan           fu = zero;
437e77315bcSPatrick Sanan         }
438e77315bcSPatrick Sanan 
439e77315bcSPatrick Sanan         f[k][j][i] = -hyhzdhx*(fe-fw) - hzhxdhy*(fn-fs) - hxhydhz*(fu-fd);
440e77315bcSPatrick Sanan       }
441e77315bcSPatrick Sanan     }
442e77315bcSPatrick Sanan   }
443*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,localX,&x));
444*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,F,&f));
445*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localX));
446*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops((22.0 + 4.0*POWFLOP)*ym*xm));
447e77315bcSPatrick Sanan   PetscFunctionReturn(0);
448e77315bcSPatrick Sanan }
449e77315bcSPatrick Sanan /* --------------------  Evaluate Jacobian F(x) --------------------- */
450e77315bcSPatrick Sanan PetscErrorCode FormJacobian(SNES snes,Vec X,Mat J,Mat jac,void *ptr)
451e77315bcSPatrick Sanan {
452e77315bcSPatrick Sanan   AppCtx         *user = (AppCtx*)ptr;
453e77315bcSPatrick Sanan   PetscInt       i,j,k,mx,my,mz,xs,ys,zs,xm,ym,zm;
454e77315bcSPatrick Sanan   PetscScalar    one = 1.0;
455e77315bcSPatrick Sanan   PetscScalar    hx,hy,hz,hxhydhz,hyhzdhx,hzhxdhy;
456e77315bcSPatrick Sanan   PetscScalar    t0,tn,ts,te,tw,an,as,ae,aw,dn,ds,de,dw;
457e77315bcSPatrick Sanan   PetscScalar    tleft,tright,beta,td,ad,dd,tu,au,du,v[7],bm1,coef;
458e77315bcSPatrick Sanan   PetscScalar    ***x,bn,bs,be,bw,bu,bd,gn,gs,ge,gw,gu,gd;
459e77315bcSPatrick Sanan   Vec            localX;
460e77315bcSPatrick Sanan   MatStencil     c[7],row;
461e77315bcSPatrick Sanan   DM             da;
462e77315bcSPatrick Sanan 
463e77315bcSPatrick Sanan   PetscFunctionBeginUser;
464*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
465*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localX));
466*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&mx,&my,&mz,0,0,0,0,0,0,0,0,0));
467e77315bcSPatrick Sanan   hx      = one/(PetscReal)(mx-1);  hy    = one/(PetscReal)(my-1);  hz = one/(PetscReal)(mz-1);
468e77315bcSPatrick Sanan   hxhydhz = hx*hy/hz;   hyhzdhx = hy*hz/hx;   hzhxdhy = hz*hx/hy;
469e77315bcSPatrick Sanan   tleft   = user->tleft;         tright = user->tright;
470e77315bcSPatrick Sanan   beta    = user->beta;          bm1    = user->bm1;              coef = user->coef;
471e77315bcSPatrick Sanan 
472e77315bcSPatrick Sanan   /* Get ghost points */
473*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
474*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
475*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm));
476*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,localX,&x));
477e77315bcSPatrick Sanan 
478e77315bcSPatrick Sanan   /* Evaluate Jacobian of function */
479e77315bcSPatrick Sanan   for (k=zs; k<zs+zm; k++) {
480e77315bcSPatrick Sanan     for (j=ys; j<ys+ym; j++) {
481e77315bcSPatrick Sanan       for (i=xs; i<xs+xm; i++) {
482e77315bcSPatrick Sanan         t0    = x[k][j][i];
483e77315bcSPatrick Sanan         row.k = k; row.j = j; row.i = i;
484e77315bcSPatrick Sanan         if (i > 0 && i < mx-1 && j > 0 && j < my-1 && k > 0 && k < mz-1) {
485e77315bcSPatrick Sanan 
486e77315bcSPatrick Sanan           /* general interior volume */
487e77315bcSPatrick Sanan 
488e77315bcSPatrick Sanan           tw = x[k][j][i-1];
489e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
490e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
491e77315bcSPatrick Sanan           /* dw = bw * aw */
492e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
493e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
494e77315bcSPatrick Sanan 
495e77315bcSPatrick Sanan           te = x[k][j][i+1];
496e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
497e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
498e77315bcSPatrick Sanan           /* de = be * ae; */
499e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
500e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
501e77315bcSPatrick Sanan 
502e77315bcSPatrick Sanan           ts = x[k][j-1][i];
503e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
504e77315bcSPatrick Sanan           bs = PetscPowScalar(as,bm1);
505e77315bcSPatrick Sanan           /* ds = bs * as; */
506e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
507e77315bcSPatrick Sanan           gs = coef*bs*(t0 - ts);
508e77315bcSPatrick Sanan 
509e77315bcSPatrick Sanan           tn = x[k][j+1][i];
510e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
511e77315bcSPatrick Sanan           bn = PetscPowScalar(an,bm1);
512e77315bcSPatrick Sanan           /* dn = bn * an; */
513e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
514e77315bcSPatrick Sanan           gn = coef*bn*(tn - t0);
515e77315bcSPatrick Sanan 
516e77315bcSPatrick Sanan           td = x[k-1][j][i];
517e77315bcSPatrick Sanan           ad = 0.5*(t0 + td);
518e77315bcSPatrick Sanan           bd = PetscPowScalar(ad,bm1);
519e77315bcSPatrick Sanan           /* dd = bd * ad; */
520e77315bcSPatrick Sanan           dd = PetscPowScalar(ad,beta);
521e77315bcSPatrick Sanan           gd = coef*bd*(t0 - td);
522e77315bcSPatrick Sanan 
523e77315bcSPatrick Sanan           tu = x[k+1][j][i];
524e77315bcSPatrick Sanan           au = 0.5*(t0 + tu);
525e77315bcSPatrick Sanan           bu = PetscPowScalar(au,bm1);
526e77315bcSPatrick Sanan           /* du = bu * au; */
527e77315bcSPatrick Sanan           du = PetscPowScalar(au,beta);
528e77315bcSPatrick Sanan           gu = coef*bu*(tu - t0);
529e77315bcSPatrick Sanan 
530e77315bcSPatrick Sanan           c[0].k = k-1; c[0].j = j; c[0].i = i; v[0]   = -hxhydhz*(dd - gd);
531e77315bcSPatrick Sanan           c[1].k = k; c[1].j = j-1; c[1].i = i;
532e77315bcSPatrick Sanan           v[1]   = -hzhxdhy*(ds - gs);
533e77315bcSPatrick Sanan           c[2].k = k; c[2].j = j; c[2].i = i-1;
534e77315bcSPatrick Sanan           v[2]   = -hyhzdhx*(dw - gw);
535e77315bcSPatrick Sanan           c[3].k = k; c[3].j = j; c[3].i = i;
536e77315bcSPatrick Sanan           v[3]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
537e77315bcSPatrick Sanan           c[4].k = k; c[4].j = j; c[4].i = i+1;
538e77315bcSPatrick Sanan           v[4]   = -hyhzdhx*(de + ge);
539e77315bcSPatrick Sanan           c[5].k = k; c[5].j = j+1; c[5].i = i;
540e77315bcSPatrick Sanan           v[5]   = -hzhxdhy*(dn + gn);
541e77315bcSPatrick Sanan           c[6].k = k+1; c[6].j = j; c[6].i = i;
542e77315bcSPatrick Sanan           v[6]   = -hxhydhz*(du + gu);
543*5f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(jac,1,&row,7,c,v,INSERT_VALUES));
544e77315bcSPatrick Sanan 
545e77315bcSPatrick Sanan         } else if (i == 0) {
546e77315bcSPatrick Sanan 
547e77315bcSPatrick Sanan           /* left-hand plane boundary */
548e77315bcSPatrick Sanan           tw = tleft;
549e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
550e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
551e77315bcSPatrick Sanan           /* dw = bw * aw */
552e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
553e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
554e77315bcSPatrick Sanan 
555e77315bcSPatrick Sanan           te = x[k][j][i+1];
556e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
557e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
558e77315bcSPatrick Sanan           /* de = be * ae; */
559e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
560e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
561e77315bcSPatrick Sanan 
562e77315bcSPatrick Sanan           /* left-hand bottom edge */
563e77315bcSPatrick Sanan           if (j == 0) {
564e77315bcSPatrick Sanan 
565e77315bcSPatrick Sanan             tn = x[k][j+1][i];
566e77315bcSPatrick Sanan             an = 0.5*(t0 + tn);
567e77315bcSPatrick Sanan             bn = PetscPowScalar(an,bm1);
568e77315bcSPatrick Sanan             /* dn = bn * an; */
569e77315bcSPatrick Sanan             dn = PetscPowScalar(an,beta);
570e77315bcSPatrick Sanan             gn = coef*bn*(tn - t0);
571e77315bcSPatrick Sanan 
572e77315bcSPatrick Sanan             /* left-hand bottom down corner */
573e77315bcSPatrick Sanan             if (k == 0) {
574e77315bcSPatrick Sanan 
575e77315bcSPatrick Sanan               tu = x[k+1][j][i];
576e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
577e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
578e77315bcSPatrick Sanan               /* du = bu * au; */
579e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
580e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
581e77315bcSPatrick Sanan 
582e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j; c[0].i = i;
583e77315bcSPatrick Sanan               v[0]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
584e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i+1;
585e77315bcSPatrick Sanan               v[1]   = -hyhzdhx*(de + ge);
586e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j+1; c[2].i = i;
587e77315bcSPatrick Sanan               v[2]   = -hzhxdhy*(dn + gn);
588e77315bcSPatrick Sanan               c[3].k = k+1; c[3].j = j; c[3].i = i;
589e77315bcSPatrick Sanan               v[3]   = -hxhydhz*(du + gu);
590*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
591e77315bcSPatrick Sanan 
592e77315bcSPatrick Sanan               /* left-hand bottom interior edge */
593e77315bcSPatrick Sanan             } else if (k < mz-1) {
594e77315bcSPatrick Sanan 
595e77315bcSPatrick Sanan               tu = x[k+1][j][i];
596e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
597e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
598e77315bcSPatrick Sanan               /* du = bu * au; */
599e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
600e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
601e77315bcSPatrick Sanan 
602e77315bcSPatrick Sanan               td = x[k-1][j][i];
603e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
604e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
605e77315bcSPatrick Sanan               /* dd = bd * ad; */
606e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
607e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
608e77315bcSPatrick Sanan 
609e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
610e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
611e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i;
612e77315bcSPatrick Sanan               v[1]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
613e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i+1;
614e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(de + ge);
615e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
616e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
617e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
618e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
619*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
620e77315bcSPatrick Sanan 
621e77315bcSPatrick Sanan               /* left-hand bottom up corner */
622e77315bcSPatrick Sanan             } else {
623e77315bcSPatrick Sanan 
624e77315bcSPatrick Sanan               td = x[k-1][j][i];
625e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
626e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
627e77315bcSPatrick Sanan               /* dd = bd * ad; */
628e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
629e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
630e77315bcSPatrick Sanan 
631e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
632e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
633e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i;
634e77315bcSPatrick Sanan               v[1]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
635e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i+1;
636e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(de + ge);
637e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
638e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
639*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
640e77315bcSPatrick Sanan             }
641e77315bcSPatrick Sanan 
642e77315bcSPatrick Sanan             /* left-hand top edge */
643e77315bcSPatrick Sanan           } else if (j == my-1) {
644e77315bcSPatrick Sanan 
645e77315bcSPatrick Sanan             ts = x[k][j-1][i];
646e77315bcSPatrick Sanan             as = 0.5*(t0 + ts);
647e77315bcSPatrick Sanan             bs = PetscPowScalar(as,bm1);
648e77315bcSPatrick Sanan             /* ds = bs * as; */
649e77315bcSPatrick Sanan             ds = PetscPowScalar(as,beta);
650e77315bcSPatrick Sanan             gs = coef*bs*(ts - t0);
651e77315bcSPatrick Sanan 
652e77315bcSPatrick Sanan             /* left-hand top down corner */
653e77315bcSPatrick Sanan             if (k == 0) {
654e77315bcSPatrick Sanan 
655e77315bcSPatrick Sanan               tu = x[k+1][j][i];
656e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
657e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
658e77315bcSPatrick Sanan               /* du = bu * au; */
659e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
660e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
661e77315bcSPatrick Sanan 
662e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j-1; c[0].i = i;
663e77315bcSPatrick Sanan               v[0]   = -hzhxdhy*(ds - gs);
664e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i;
665e77315bcSPatrick Sanan               v[1]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
666e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i+1;
667e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(de + ge);
668e77315bcSPatrick Sanan               c[3].k = k+1; c[3].j = j; c[3].i = i;
669e77315bcSPatrick Sanan               v[3]   = -hxhydhz*(du + gu);
670*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
671e77315bcSPatrick Sanan 
672e77315bcSPatrick Sanan               /* left-hand top interior edge */
673e77315bcSPatrick Sanan             } else if (k < mz-1) {
674e77315bcSPatrick Sanan 
675e77315bcSPatrick Sanan               tu = x[k+1][j][i];
676e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
677e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
678e77315bcSPatrick Sanan               /* du = bu * au; */
679e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
680e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
681e77315bcSPatrick Sanan 
682e77315bcSPatrick Sanan               td = x[k-1][j][i];
683e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
684e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
685e77315bcSPatrick Sanan               /* dd = bd * ad; */
686e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
687e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
688e77315bcSPatrick Sanan 
689e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
690e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
691e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
692e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
693e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
694e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
695e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i+1;
696e77315bcSPatrick Sanan               v[3]   = -hyhzdhx*(de + ge);
697e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
698e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
699*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
700e77315bcSPatrick Sanan 
701e77315bcSPatrick Sanan               /* left-hand top up corner */
702e77315bcSPatrick Sanan             } else {
703e77315bcSPatrick Sanan 
704e77315bcSPatrick Sanan               td = x[k-1][j][i];
705e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
706e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
707e77315bcSPatrick Sanan               /* dd = bd * ad; */
708e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
709e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
710e77315bcSPatrick Sanan 
711e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
712e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
713e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
714e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
715e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
716e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
717e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i+1;
718e77315bcSPatrick Sanan               v[3]   = -hyhzdhx*(de + ge);
719*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
720e77315bcSPatrick Sanan             }
721e77315bcSPatrick Sanan 
722e77315bcSPatrick Sanan           } else {
723e77315bcSPatrick Sanan 
724e77315bcSPatrick Sanan             ts = x[k][j-1][i];
725e77315bcSPatrick Sanan             as = 0.5*(t0 + ts);
726e77315bcSPatrick Sanan             bs = PetscPowScalar(as,bm1);
727e77315bcSPatrick Sanan             /* ds = bs * as; */
728e77315bcSPatrick Sanan             ds = PetscPowScalar(as,beta);
729e77315bcSPatrick Sanan             gs = coef*bs*(t0 - ts);
730e77315bcSPatrick Sanan 
731e77315bcSPatrick Sanan             tn = x[k][j+1][i];
732e77315bcSPatrick Sanan             an = 0.5*(t0 + tn);
733e77315bcSPatrick Sanan             bn = PetscPowScalar(an,bm1);
734e77315bcSPatrick Sanan             /* dn = bn * an; */
735e77315bcSPatrick Sanan             dn = PetscPowScalar(an,beta);
736e77315bcSPatrick Sanan             gn = coef*bn*(tn - t0);
737e77315bcSPatrick Sanan 
738e77315bcSPatrick Sanan             /* left-hand down interior edge */
739e77315bcSPatrick Sanan             if (k == 0) {
740e77315bcSPatrick Sanan 
741e77315bcSPatrick Sanan               tu = x[k+1][j][i];
742e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
743e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
744e77315bcSPatrick Sanan               /* du = bu * au; */
745e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
746e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
747e77315bcSPatrick Sanan 
748e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j-1; c[0].i = i;
749e77315bcSPatrick Sanan               v[0]   = -hzhxdhy*(ds - gs);
750e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i;
751e77315bcSPatrick Sanan               v[1]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
752e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i+1;
753e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(de + ge);
754e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
755e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
756e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
757e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
758*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
759e77315bcSPatrick Sanan 
760e77315bcSPatrick Sanan             } else if (k == mz-1) { /* left-hand up interior edge */
761e77315bcSPatrick Sanan 
762e77315bcSPatrick Sanan               td = x[k-1][j][i];
763e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
764e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
765e77315bcSPatrick Sanan               /* dd = bd * ad; */
766e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
767e77315bcSPatrick Sanan               gd = coef*bd*(t0 - td);
768e77315bcSPatrick Sanan 
769e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
770e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
771e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
772e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
773e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
774e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
775e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i+1;
776e77315bcSPatrick Sanan               v[3]   = -hyhzdhx*(de + ge);
777e77315bcSPatrick Sanan               c[4].k = k; c[4].j = j+1; c[4].i = i;
778e77315bcSPatrick Sanan               v[4]   = -hzhxdhy*(dn + gn);
779*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
780e77315bcSPatrick Sanan             } else { /* left-hand interior plane */
781e77315bcSPatrick Sanan 
782e77315bcSPatrick Sanan               td = x[k-1][j][i];
783e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
784e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
785e77315bcSPatrick Sanan               /* dd = bd * ad; */
786e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
787e77315bcSPatrick Sanan               gd = coef*bd*(t0 - td);
788e77315bcSPatrick Sanan 
789e77315bcSPatrick Sanan               tu = x[k+1][j][i];
790e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
791e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
792e77315bcSPatrick Sanan               /* du = bu * au; */
793e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
794e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
795e77315bcSPatrick Sanan 
796e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
797e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
798e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
799e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
800e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
801e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
802e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i+1;
803e77315bcSPatrick Sanan               v[3]   = -hyhzdhx*(de + ge);
804e77315bcSPatrick Sanan               c[4].k = k; c[4].j = j+1; c[4].i = i;
805e77315bcSPatrick Sanan               v[4]   = -hzhxdhy*(dn + gn);
806e77315bcSPatrick Sanan               c[5].k = k+1; c[5].j = j; c[5].i = i;
807e77315bcSPatrick Sanan               v[5]   = -hxhydhz*(du + gu);
808*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
809e77315bcSPatrick Sanan             }
810e77315bcSPatrick Sanan           }
811e77315bcSPatrick Sanan 
812e77315bcSPatrick Sanan         } else if (i == mx-1) {
813e77315bcSPatrick Sanan 
814e77315bcSPatrick Sanan           /* right-hand plane boundary */
815e77315bcSPatrick Sanan           tw = x[k][j][i-1];
816e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
817e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
818e77315bcSPatrick Sanan           /* dw = bw * aw */
819e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
820e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
821e77315bcSPatrick Sanan 
822e77315bcSPatrick Sanan           te = tright;
823e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
824e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
825e77315bcSPatrick Sanan           /* de = be * ae; */
826e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
827e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
828e77315bcSPatrick Sanan 
829e77315bcSPatrick Sanan           /* right-hand bottom edge */
830e77315bcSPatrick Sanan           if (j == 0) {
831e77315bcSPatrick Sanan 
832e77315bcSPatrick Sanan             tn = x[k][j+1][i];
833e77315bcSPatrick Sanan             an = 0.5*(t0 + tn);
834e77315bcSPatrick Sanan             bn = PetscPowScalar(an,bm1);
835e77315bcSPatrick Sanan             /* dn = bn * an; */
836e77315bcSPatrick Sanan             dn = PetscPowScalar(an,beta);
837e77315bcSPatrick Sanan             gn = coef*bn*(tn - t0);
838e77315bcSPatrick Sanan 
839e77315bcSPatrick Sanan             /* right-hand bottom down corner */
840e77315bcSPatrick Sanan             if (k == 0) {
841e77315bcSPatrick Sanan 
842e77315bcSPatrick Sanan               tu = x[k+1][j][i];
843e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
844e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
845e77315bcSPatrick Sanan               /* du = bu * au; */
846e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
847e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
848e77315bcSPatrick Sanan 
849e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j; c[0].i = i-1;
850e77315bcSPatrick Sanan               v[0]   = -hyhzdhx*(dw - gw);
851e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i;
852e77315bcSPatrick Sanan               v[1]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
853e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j+1; c[2].i = i;
854e77315bcSPatrick Sanan               v[2]   = -hzhxdhy*(dn + gn);
855e77315bcSPatrick Sanan               c[3].k = k+1; c[3].j = j; c[3].i = i;
856e77315bcSPatrick Sanan               v[3]   = -hxhydhz*(du + gu);
857*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
858e77315bcSPatrick Sanan 
859e77315bcSPatrick Sanan               /* right-hand bottom interior edge */
860e77315bcSPatrick Sanan             } else if (k < mz-1) {
861e77315bcSPatrick Sanan 
862e77315bcSPatrick Sanan               tu = x[k+1][j][i];
863e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
864e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
865e77315bcSPatrick Sanan               /* du = bu * au; */
866e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
867e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
868e77315bcSPatrick Sanan 
869e77315bcSPatrick Sanan               td = x[k-1][j][i];
870e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
871e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
872e77315bcSPatrick Sanan               /* dd = bd * ad; */
873e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
874e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
875e77315bcSPatrick Sanan 
876e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
877e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
878e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i-1;
879e77315bcSPatrick Sanan               v[1]   = -hyhzdhx*(dw - gw);
880e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
881e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
882e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
883e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
884e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
885e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
886*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
887e77315bcSPatrick Sanan 
888e77315bcSPatrick Sanan               /* right-hand bottom up corner */
889e77315bcSPatrick Sanan             } else {
890e77315bcSPatrick Sanan 
891e77315bcSPatrick Sanan               td = x[k-1][j][i];
892e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
893e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
894e77315bcSPatrick Sanan               /* dd = bd * ad; */
895e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
896e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
897e77315bcSPatrick Sanan 
898e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
899e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
900e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i-1;
901e77315bcSPatrick Sanan               v[1]   = -hyhzdhx*(dw - gw);
902e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
903e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
904e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
905e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
906*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
907e77315bcSPatrick Sanan             }
908e77315bcSPatrick Sanan 
909e77315bcSPatrick Sanan             /* right-hand top edge */
910e77315bcSPatrick Sanan           } else if (j == my-1) {
911e77315bcSPatrick Sanan 
912e77315bcSPatrick Sanan             ts = x[k][j-1][i];
913e77315bcSPatrick Sanan             as = 0.5*(t0 + ts);
914e77315bcSPatrick Sanan             bs = PetscPowScalar(as,bm1);
915e77315bcSPatrick Sanan             /* ds = bs * as; */
916e77315bcSPatrick Sanan             ds = PetscPowScalar(as,beta);
917e77315bcSPatrick Sanan             gs = coef*bs*(ts - t0);
918e77315bcSPatrick Sanan 
919e77315bcSPatrick Sanan             /* right-hand top down corner */
920e77315bcSPatrick Sanan             if (k == 0) {
921e77315bcSPatrick Sanan 
922e77315bcSPatrick Sanan               tu = x[k+1][j][i];
923e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
924e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
925e77315bcSPatrick Sanan               /* du = bu * au; */
926e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
927e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
928e77315bcSPatrick Sanan 
929e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j-1; c[0].i = i;
930e77315bcSPatrick Sanan               v[0]   = -hzhxdhy*(ds - gs);
931e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i-1;
932e77315bcSPatrick Sanan               v[1]   = -hyhzdhx*(dw - gw);
933e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
934e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
935e77315bcSPatrick Sanan               c[3].k = k+1; c[3].j = j; c[3].i = i;
936e77315bcSPatrick Sanan               v[3]   = -hxhydhz*(du + gu);
937*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
938e77315bcSPatrick Sanan 
939e77315bcSPatrick Sanan               /* right-hand top interior edge */
940e77315bcSPatrick Sanan             } else if (k < mz-1) {
941e77315bcSPatrick Sanan 
942e77315bcSPatrick Sanan               tu = x[k+1][j][i];
943e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
944e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
945e77315bcSPatrick Sanan               /* du = bu * au; */
946e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
947e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
948e77315bcSPatrick Sanan 
949e77315bcSPatrick Sanan               td = x[k-1][j][i];
950e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
951e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
952e77315bcSPatrick Sanan               /* dd = bd * ad; */
953e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
954e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
955e77315bcSPatrick Sanan 
956e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
957e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
958e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
959e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
960e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i-1;
961e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(dw - gw);
962e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i;
963e77315bcSPatrick Sanan               v[3]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
964e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
965e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
966*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
967e77315bcSPatrick Sanan 
968e77315bcSPatrick Sanan               /* right-hand top up corner */
969e77315bcSPatrick Sanan             } else {
970e77315bcSPatrick Sanan 
971e77315bcSPatrick Sanan               td = x[k-1][j][i];
972e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
973e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
974e77315bcSPatrick Sanan               /* dd = bd * ad; */
975e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
976e77315bcSPatrick Sanan               gd = coef*bd*(td - t0);
977e77315bcSPatrick Sanan 
978e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
979e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
980e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
981e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
982e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i-1;
983e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(dw - gw);
984e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i;
985e77315bcSPatrick Sanan               v[3]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
986*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,4,c,v,INSERT_VALUES));
987e77315bcSPatrick Sanan             }
988e77315bcSPatrick Sanan 
989e77315bcSPatrick Sanan           } else {
990e77315bcSPatrick Sanan 
991e77315bcSPatrick Sanan             ts = x[k][j-1][i];
992e77315bcSPatrick Sanan             as = 0.5*(t0 + ts);
993e77315bcSPatrick Sanan             bs = PetscPowScalar(as,bm1);
994e77315bcSPatrick Sanan             /* ds = bs * as; */
995e77315bcSPatrick Sanan             ds = PetscPowScalar(as,beta);
996e77315bcSPatrick Sanan             gs = coef*bs*(t0 - ts);
997e77315bcSPatrick Sanan 
998e77315bcSPatrick Sanan             tn = x[k][j+1][i];
999e77315bcSPatrick Sanan             an = 0.5*(t0 + tn);
1000e77315bcSPatrick Sanan             bn = PetscPowScalar(an,bm1);
1001e77315bcSPatrick Sanan             /* dn = bn * an; */
1002e77315bcSPatrick Sanan             dn = PetscPowScalar(an,beta);
1003e77315bcSPatrick Sanan             gn = coef*bn*(tn - t0);
1004e77315bcSPatrick Sanan 
1005e77315bcSPatrick Sanan             /* right-hand down interior edge */
1006e77315bcSPatrick Sanan             if (k == 0) {
1007e77315bcSPatrick Sanan 
1008e77315bcSPatrick Sanan               tu = x[k+1][j][i];
1009e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
1010e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
1011e77315bcSPatrick Sanan               /* du = bu * au; */
1012e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
1013e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
1014e77315bcSPatrick Sanan 
1015e77315bcSPatrick Sanan               c[0].k = k; c[0].j = j-1; c[0].i = i;
1016e77315bcSPatrick Sanan               v[0]   = -hzhxdhy*(ds - gs);
1017e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j; c[1].i = i-1;
1018e77315bcSPatrick Sanan               v[1]   = -hyhzdhx*(dw - gw);
1019e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i;
1020e77315bcSPatrick Sanan               v[2]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
1021e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j+1; c[3].i = i;
1022e77315bcSPatrick Sanan               v[3]   = -hzhxdhy*(dn + gn);
1023e77315bcSPatrick Sanan               c[4].k = k+1; c[4].j = j; c[4].i = i;
1024e77315bcSPatrick Sanan               v[4]   = -hxhydhz*(du + gu);
1025*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1026e77315bcSPatrick Sanan 
1027e77315bcSPatrick Sanan             } else if (k == mz-1) { /* right-hand up interior edge */
1028e77315bcSPatrick Sanan 
1029e77315bcSPatrick Sanan               td = x[k-1][j][i];
1030e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
1031e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
1032e77315bcSPatrick Sanan               /* dd = bd * ad; */
1033e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
1034e77315bcSPatrick Sanan               gd = coef*bd*(t0 - td);
1035e77315bcSPatrick Sanan 
1036e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
1037e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
1038e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
1039e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
1040e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i-1;
1041e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(dw - gw);
1042e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i;
1043e77315bcSPatrick Sanan               v[3]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
1044e77315bcSPatrick Sanan               c[4].k = k; c[4].j = j+1; c[4].i = i;
1045e77315bcSPatrick Sanan               v[4]   = -hzhxdhy*(dn + gn);
1046*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1047e77315bcSPatrick Sanan 
1048e77315bcSPatrick Sanan             } else { /* right-hand interior plane */
1049e77315bcSPatrick Sanan 
1050e77315bcSPatrick Sanan               td = x[k-1][j][i];
1051e77315bcSPatrick Sanan               ad = 0.5*(t0 + td);
1052e77315bcSPatrick Sanan               bd = PetscPowScalar(ad,bm1);
1053e77315bcSPatrick Sanan               /* dd = bd * ad; */
1054e77315bcSPatrick Sanan               dd = PetscPowScalar(ad,beta);
1055e77315bcSPatrick Sanan               gd = coef*bd*(t0 - td);
1056e77315bcSPatrick Sanan 
1057e77315bcSPatrick Sanan               tu = x[k+1][j][i];
1058e77315bcSPatrick Sanan               au = 0.5*(t0 + tu);
1059e77315bcSPatrick Sanan               bu = PetscPowScalar(au,bm1);
1060e77315bcSPatrick Sanan               /* du = bu * au; */
1061e77315bcSPatrick Sanan               du = PetscPowScalar(au,beta);
1062e77315bcSPatrick Sanan               gu = coef*bu*(tu - t0);
1063e77315bcSPatrick Sanan 
1064e77315bcSPatrick Sanan               c[0].k = k-1; c[0].j = j; c[0].i = i;
1065e77315bcSPatrick Sanan               v[0]   = -hxhydhz*(dd - gd);
1066e77315bcSPatrick Sanan               c[1].k = k; c[1].j = j-1; c[1].i = i;
1067e77315bcSPatrick Sanan               v[1]   = -hzhxdhy*(ds - gs);
1068e77315bcSPatrick Sanan               c[2].k = k; c[2].j = j; c[2].i = i-1;
1069e77315bcSPatrick Sanan               v[2]   = -hyhzdhx*(dw - gw);
1070e77315bcSPatrick Sanan               c[3].k = k; c[3].j = j; c[3].i = i;
1071e77315bcSPatrick Sanan               v[3]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
1072e77315bcSPatrick Sanan               c[4].k = k; c[4].j = j+1; c[4].i = i;
1073e77315bcSPatrick Sanan               v[4]   = -hzhxdhy*(dn + gn);
1074e77315bcSPatrick Sanan               c[5].k = k+1; c[5].j = j; c[5].i = i;
1075e77315bcSPatrick Sanan               v[5]   = -hxhydhz*(du + gu);
1076*5f80ce2aSJacob Faibussowitsch               CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
1077e77315bcSPatrick Sanan             }
1078e77315bcSPatrick Sanan           }
1079e77315bcSPatrick Sanan 
1080e77315bcSPatrick Sanan         } else if (j == 0) {
1081e77315bcSPatrick Sanan 
1082e77315bcSPatrick Sanan           tw = x[k][j][i-1];
1083e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
1084e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
1085e77315bcSPatrick Sanan           /* dw = bw * aw */
1086e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
1087e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
1088e77315bcSPatrick Sanan 
1089e77315bcSPatrick Sanan           te = x[k][j][i+1];
1090e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
1091e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
1092e77315bcSPatrick Sanan           /* de = be * ae; */
1093e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
1094e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
1095e77315bcSPatrick Sanan 
1096e77315bcSPatrick Sanan           tn = x[k][j+1][i];
1097e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
1098e77315bcSPatrick Sanan           bn = PetscPowScalar(an,bm1);
1099e77315bcSPatrick Sanan           /* dn = bn * an; */
1100e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
1101e77315bcSPatrick Sanan           gn = coef*bn*(tn - t0);
1102e77315bcSPatrick Sanan 
1103e77315bcSPatrick Sanan           /* bottom down interior edge */
1104e77315bcSPatrick Sanan           if (k == 0) {
1105e77315bcSPatrick Sanan 
1106e77315bcSPatrick Sanan             tu = x[k+1][j][i];
1107e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
1108e77315bcSPatrick Sanan             bu = PetscPowScalar(au,bm1);
1109e77315bcSPatrick Sanan             /* du = bu * au; */
1110e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
1111e77315bcSPatrick Sanan             gu = coef*bu*(tu - t0);
1112e77315bcSPatrick Sanan 
1113e77315bcSPatrick Sanan             c[0].k = k; c[0].j = j; c[0].i = i-1;
1114e77315bcSPatrick Sanan             v[0]   = -hyhzdhx*(dw - gw);
1115e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j; c[1].i = i;
1116e77315bcSPatrick Sanan             v[1]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
1117e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i+1;
1118e77315bcSPatrick Sanan             v[2]   = -hyhzdhx*(de + ge);
1119e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j+1; c[3].i = i;
1120e77315bcSPatrick Sanan             v[3]   = -hzhxdhy*(dn + gn);
1121e77315bcSPatrick Sanan             c[4].k = k+1; c[4].j = j; c[4].i = i;
1122e77315bcSPatrick Sanan             v[4]   = -hxhydhz*(du + gu);
1123*5f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1124e77315bcSPatrick Sanan 
1125e77315bcSPatrick Sanan           } else if (k == mz-1) { /* bottom up interior edge */
1126e77315bcSPatrick Sanan 
1127e77315bcSPatrick Sanan             td = x[k-1][j][i];
1128e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
1129e77315bcSPatrick Sanan             bd = PetscPowScalar(ad,bm1);
1130e77315bcSPatrick Sanan             /* dd = bd * ad; */
1131e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
1132e77315bcSPatrick Sanan             gd = coef*bd*(td - t0);
1133e77315bcSPatrick Sanan 
1134e77315bcSPatrick Sanan             c[0].k = k-1; c[0].j = j; c[0].i = i;
1135e77315bcSPatrick Sanan             v[0]   = -hxhydhz*(dd - gd);
1136e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j; c[1].i = i-1;
1137e77315bcSPatrick Sanan             v[1]   = -hyhzdhx*(dw - gw);
1138e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i;
1139e77315bcSPatrick Sanan             v[2]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
1140e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j; c[3].i = i+1;
1141e77315bcSPatrick Sanan             v[3]   = -hyhzdhx*(de + ge);
1142e77315bcSPatrick Sanan             c[4].k = k; c[4].j = j+1; c[4].i = i;
1143e77315bcSPatrick Sanan             v[4]   = -hzhxdhy*(dn + gn);
1144*5f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1145e77315bcSPatrick Sanan 
1146e77315bcSPatrick Sanan           } else { /* bottom interior plane */
1147e77315bcSPatrick Sanan 
1148e77315bcSPatrick Sanan             tu = x[k+1][j][i];
1149e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
1150e77315bcSPatrick Sanan             bu = PetscPowScalar(au,bm1);
1151e77315bcSPatrick Sanan             /* du = bu * au; */
1152e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
1153e77315bcSPatrick Sanan             gu = coef*bu*(tu - t0);
1154e77315bcSPatrick Sanan 
1155e77315bcSPatrick Sanan             td = x[k-1][j][i];
1156e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
1157e77315bcSPatrick Sanan             bd = PetscPowScalar(ad,bm1);
1158e77315bcSPatrick Sanan             /* dd = bd * ad; */
1159e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
1160e77315bcSPatrick Sanan             gd = coef*bd*(td - t0);
1161e77315bcSPatrick Sanan 
1162e77315bcSPatrick Sanan             c[0].k = k-1; c[0].j = j; c[0].i = i;
1163e77315bcSPatrick Sanan             v[0]   = -hxhydhz*(dd - gd);
1164e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j; c[1].i = i-1;
1165e77315bcSPatrick Sanan             v[1]   = -hyhzdhx*(dw - gw);
1166e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i;
1167e77315bcSPatrick Sanan             v[2]   =   hzhxdhy*(dn - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
1168e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j; c[3].i = i+1;
1169e77315bcSPatrick Sanan             v[3]   = -hyhzdhx*(de + ge);
1170e77315bcSPatrick Sanan             c[4].k = k; c[4].j = j+1; c[4].i = i;
1171e77315bcSPatrick Sanan             v[4]   = -hzhxdhy*(dn + gn);
1172e77315bcSPatrick Sanan             c[5].k = k+1; c[5].j = j; c[5].i = i;
1173e77315bcSPatrick Sanan             v[5]   = -hxhydhz*(du + gu);
1174*5f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
1175e77315bcSPatrick Sanan           }
1176e77315bcSPatrick Sanan 
1177e77315bcSPatrick Sanan         } else if (j == my-1) {
1178e77315bcSPatrick Sanan 
1179e77315bcSPatrick Sanan           tw = x[k][j][i-1];
1180e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
1181e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
1182e77315bcSPatrick Sanan           /* dw = bw * aw */
1183e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
1184e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
1185e77315bcSPatrick Sanan 
1186e77315bcSPatrick Sanan           te = x[k][j][i+1];
1187e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
1188e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
1189e77315bcSPatrick Sanan           /* de = be * ae; */
1190e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
1191e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
1192e77315bcSPatrick Sanan 
1193e77315bcSPatrick Sanan           ts = x[k][j-1][i];
1194e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
1195e77315bcSPatrick Sanan           bs = PetscPowScalar(as,bm1);
1196e77315bcSPatrick Sanan           /* ds = bs * as; */
1197e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
1198e77315bcSPatrick Sanan           gs = coef*bs*(t0 - ts);
1199e77315bcSPatrick Sanan 
1200e77315bcSPatrick Sanan           /* top down interior edge */
1201e77315bcSPatrick Sanan           if (k == 0) {
1202e77315bcSPatrick Sanan 
1203e77315bcSPatrick Sanan             tu = x[k+1][j][i];
1204e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
1205e77315bcSPatrick Sanan             bu = PetscPowScalar(au,bm1);
1206e77315bcSPatrick Sanan             /* du = bu * au; */
1207e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
1208e77315bcSPatrick Sanan             gu = coef*bu*(tu - t0);
1209e77315bcSPatrick Sanan 
1210e77315bcSPatrick Sanan             c[0].k = k; c[0].j = j-1; c[0].i = i;
1211e77315bcSPatrick Sanan             v[0]   = -hzhxdhy*(ds - gs);
1212e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j; c[1].i = i-1;
1213e77315bcSPatrick Sanan             v[1]   = -hyhzdhx*(dw - gw);
1214e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i;
1215e77315bcSPatrick Sanan             v[2]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
1216e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j; c[3].i = i+1;
1217e77315bcSPatrick Sanan             v[3]   = -hyhzdhx*(de + ge);
1218e77315bcSPatrick Sanan             c[4].k = k+1; c[4].j = j; c[4].i = i;
1219e77315bcSPatrick Sanan             v[4]   = -hxhydhz*(du + gu);
1220*5f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1221e77315bcSPatrick Sanan 
1222e77315bcSPatrick Sanan           } else if (k == mz-1) { /* top up interior edge */
1223e77315bcSPatrick Sanan 
1224e77315bcSPatrick Sanan             td = x[k-1][j][i];
1225e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
1226e77315bcSPatrick Sanan             bd = PetscPowScalar(ad,bm1);
1227e77315bcSPatrick Sanan             /* dd = bd * ad; */
1228e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
1229e77315bcSPatrick Sanan             gd = coef*bd*(td - t0);
1230e77315bcSPatrick Sanan 
1231e77315bcSPatrick Sanan             c[0].k = k-1; c[0].j = j; c[0].i = i;
1232e77315bcSPatrick Sanan             v[0]   = -hxhydhz*(dd - gd);
1233e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j-1; c[1].i = i;
1234e77315bcSPatrick Sanan             v[1]   = -hzhxdhy*(ds - gs);
1235e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i-1;
1236e77315bcSPatrick Sanan             v[2]   = -hyhzdhx*(dw - gw);
1237e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j; c[3].i = i;
1238e77315bcSPatrick Sanan             v[3]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
1239e77315bcSPatrick Sanan             c[4].k = k; c[4].j = j; c[4].i = i+1;
1240e77315bcSPatrick Sanan             v[4]   = -hyhzdhx*(de + ge);
1241*5f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,5,c,v,INSERT_VALUES));
1242e77315bcSPatrick Sanan 
1243e77315bcSPatrick Sanan           } else { /* top interior plane */
1244e77315bcSPatrick Sanan 
1245e77315bcSPatrick Sanan             tu = x[k+1][j][i];
1246e77315bcSPatrick Sanan             au = 0.5*(t0 + tu);
1247e77315bcSPatrick Sanan             bu = PetscPowScalar(au,bm1);
1248e77315bcSPatrick Sanan             /* du = bu * au; */
1249e77315bcSPatrick Sanan             du = PetscPowScalar(au,beta);
1250e77315bcSPatrick Sanan             gu = coef*bu*(tu - t0);
1251e77315bcSPatrick Sanan 
1252e77315bcSPatrick Sanan             td = x[k-1][j][i];
1253e77315bcSPatrick Sanan             ad = 0.5*(t0 + td);
1254e77315bcSPatrick Sanan             bd = PetscPowScalar(ad,bm1);
1255e77315bcSPatrick Sanan             /* dd = bd * ad; */
1256e77315bcSPatrick Sanan             dd = PetscPowScalar(ad,beta);
1257e77315bcSPatrick Sanan             gd = coef*bd*(td - t0);
1258e77315bcSPatrick Sanan 
1259e77315bcSPatrick Sanan             c[0].k = k-1; c[0].j = j; c[0].i = i;
1260e77315bcSPatrick Sanan             v[0]   = -hxhydhz*(dd - gd);
1261e77315bcSPatrick Sanan             c[1].k = k; c[1].j = j-1; c[1].i = i;
1262e77315bcSPatrick Sanan             v[1]   = -hzhxdhy*(ds - gs);
1263e77315bcSPatrick Sanan             c[2].k = k; c[2].j = j; c[2].i = i-1;
1264e77315bcSPatrick Sanan             v[2]   = -hyhzdhx*(dw - gw);
1265e77315bcSPatrick Sanan             c[3].k = k; c[3].j = j; c[3].i = i;
1266e77315bcSPatrick Sanan             v[3]   =   hzhxdhy*(ds + gs) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + du + gd - gu);
1267e77315bcSPatrick Sanan             c[4].k = k; c[4].j = j; c[4].i = i+1;
1268e77315bcSPatrick Sanan             v[4]   = -hyhzdhx*(de + ge);
1269e77315bcSPatrick Sanan             c[5].k = k+1; c[5].j = j; c[5].i = i;
1270e77315bcSPatrick Sanan             v[5]   = -hxhydhz*(du + gu);
1271*5f80ce2aSJacob Faibussowitsch             CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
1272e77315bcSPatrick Sanan           }
1273e77315bcSPatrick Sanan 
1274e77315bcSPatrick Sanan         } else if (k == 0) {
1275e77315bcSPatrick Sanan 
1276e77315bcSPatrick Sanan           /* down interior plane */
1277e77315bcSPatrick Sanan 
1278e77315bcSPatrick Sanan           tw = x[k][j][i-1];
1279e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
1280e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
1281e77315bcSPatrick Sanan           /* dw = bw * aw */
1282e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
1283e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
1284e77315bcSPatrick Sanan 
1285e77315bcSPatrick Sanan           te = x[k][j][i+1];
1286e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
1287e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
1288e77315bcSPatrick Sanan           /* de = be * ae; */
1289e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
1290e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
1291e77315bcSPatrick Sanan 
1292e77315bcSPatrick Sanan           ts = x[k][j-1][i];
1293e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
1294e77315bcSPatrick Sanan           bs = PetscPowScalar(as,bm1);
1295e77315bcSPatrick Sanan           /* ds = bs * as; */
1296e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
1297e77315bcSPatrick Sanan           gs = coef*bs*(t0 - ts);
1298e77315bcSPatrick Sanan 
1299e77315bcSPatrick Sanan           tn = x[k][j+1][i];
1300e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
1301e77315bcSPatrick Sanan           bn = PetscPowScalar(an,bm1);
1302e77315bcSPatrick Sanan           /* dn = bn * an; */
1303e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
1304e77315bcSPatrick Sanan           gn = coef*bn*(tn - t0);
1305e77315bcSPatrick Sanan 
1306e77315bcSPatrick Sanan           tu = x[k+1][j][i];
1307e77315bcSPatrick Sanan           au = 0.5*(t0 + tu);
1308e77315bcSPatrick Sanan           bu = PetscPowScalar(au,bm1);
1309e77315bcSPatrick Sanan           /* du = bu * au; */
1310e77315bcSPatrick Sanan           du = PetscPowScalar(au,beta);
1311e77315bcSPatrick Sanan           gu = coef*bu*(tu - t0);
1312e77315bcSPatrick Sanan 
1313e77315bcSPatrick Sanan           c[0].k = k; c[0].j = j-1; c[0].i = i;
1314e77315bcSPatrick Sanan           v[0]   = -hzhxdhy*(ds - gs);
1315e77315bcSPatrick Sanan           c[1].k = k; c[1].j = j; c[1].i = i-1;
1316e77315bcSPatrick Sanan           v[1]   = -hyhzdhx*(dw - gw);
1317e77315bcSPatrick Sanan           c[2].k = k; c[2].j = j; c[2].i = i;
1318e77315bcSPatrick Sanan           v[2]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(du - gu);
1319e77315bcSPatrick Sanan           c[3].k = k; c[3].j = j; c[3].i = i+1;
1320e77315bcSPatrick Sanan           v[3]   = -hyhzdhx*(de + ge);
1321e77315bcSPatrick Sanan           c[4].k = k; c[4].j = j+1; c[4].i = i;
1322e77315bcSPatrick Sanan           v[4]   = -hzhxdhy*(dn + gn);
1323e77315bcSPatrick Sanan           c[5].k = k+1; c[5].j = j; c[5].i = i;
1324e77315bcSPatrick Sanan           v[5]   = -hxhydhz*(du + gu);
1325*5f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
1326e77315bcSPatrick Sanan 
1327e77315bcSPatrick Sanan         } else if (k == mz-1) {
1328e77315bcSPatrick Sanan 
1329e77315bcSPatrick Sanan           /* up interior plane */
1330e77315bcSPatrick Sanan 
1331e77315bcSPatrick Sanan           tw = x[k][j][i-1];
1332e77315bcSPatrick Sanan           aw = 0.5*(t0 + tw);
1333e77315bcSPatrick Sanan           bw = PetscPowScalar(aw,bm1);
1334e77315bcSPatrick Sanan           /* dw = bw * aw */
1335e77315bcSPatrick Sanan           dw = PetscPowScalar(aw,beta);
1336e77315bcSPatrick Sanan           gw = coef*bw*(t0 - tw);
1337e77315bcSPatrick Sanan 
1338e77315bcSPatrick Sanan           te = x[k][j][i+1];
1339e77315bcSPatrick Sanan           ae = 0.5*(t0 + te);
1340e77315bcSPatrick Sanan           be = PetscPowScalar(ae,bm1);
1341e77315bcSPatrick Sanan           /* de = be * ae; */
1342e77315bcSPatrick Sanan           de = PetscPowScalar(ae,beta);
1343e77315bcSPatrick Sanan           ge = coef*be*(te - t0);
1344e77315bcSPatrick Sanan 
1345e77315bcSPatrick Sanan           ts = x[k][j-1][i];
1346e77315bcSPatrick Sanan           as = 0.5*(t0 + ts);
1347e77315bcSPatrick Sanan           bs = PetscPowScalar(as,bm1);
1348e77315bcSPatrick Sanan           /* ds = bs * as; */
1349e77315bcSPatrick Sanan           ds = PetscPowScalar(as,beta);
1350e77315bcSPatrick Sanan           gs = coef*bs*(t0 - ts);
1351e77315bcSPatrick Sanan 
1352e77315bcSPatrick Sanan           tn = x[k][j+1][i];
1353e77315bcSPatrick Sanan           an = 0.5*(t0 + tn);
1354e77315bcSPatrick Sanan           bn = PetscPowScalar(an,bm1);
1355e77315bcSPatrick Sanan           /* dn = bn * an; */
1356e77315bcSPatrick Sanan           dn = PetscPowScalar(an,beta);
1357e77315bcSPatrick Sanan           gn = coef*bn*(tn - t0);
1358e77315bcSPatrick Sanan 
1359e77315bcSPatrick Sanan           td = x[k-1][j][i];
1360e77315bcSPatrick Sanan           ad = 0.5*(t0 + td);
1361e77315bcSPatrick Sanan           bd = PetscPowScalar(ad,bm1);
1362e77315bcSPatrick Sanan           /* dd = bd * ad; */
1363e77315bcSPatrick Sanan           dd = PetscPowScalar(ad,beta);
1364e77315bcSPatrick Sanan           gd = coef*bd*(t0 - td);
1365e77315bcSPatrick Sanan 
1366e77315bcSPatrick Sanan           c[0].k = k-1; c[0].j = j; c[0].i = i;
1367e77315bcSPatrick Sanan           v[0]   = -hxhydhz*(dd - gd);
1368e77315bcSPatrick Sanan           c[1].k = k; c[1].j = j-1; c[1].i = i;
1369e77315bcSPatrick Sanan           v[1]   = -hzhxdhy*(ds - gs);
1370e77315bcSPatrick Sanan           c[2].k = k; c[2].j = j; c[2].i = i-1;
1371e77315bcSPatrick Sanan           v[2]   = -hyhzdhx*(dw - gw);
1372e77315bcSPatrick Sanan           c[3].k = k; c[3].j = j; c[3].i = i;
1373e77315bcSPatrick Sanan           v[3]   =   hzhxdhy*(ds + dn + gs - gn) + hyhzdhx*(dw + de + gw - ge) + hxhydhz*(dd + gd);
1374e77315bcSPatrick Sanan           c[4].k = k; c[4].j = j; c[4].i = i+1;
1375e77315bcSPatrick Sanan           v[4]   = -hyhzdhx*(de + ge);
1376e77315bcSPatrick Sanan           c[5].k = k; c[5].j = j+1; c[5].i = i;
1377e77315bcSPatrick Sanan           v[5]   = -hzhxdhy*(dn + gn);
1378*5f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(jac,1,&row,6,c,v,INSERT_VALUES));
1379e77315bcSPatrick Sanan         }
1380e77315bcSPatrick Sanan       }
1381e77315bcSPatrick Sanan     }
1382e77315bcSPatrick Sanan   }
1383*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
1384*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
1385e77315bcSPatrick Sanan   if (jac != J) {
1386*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
1387*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
1388e77315bcSPatrick Sanan   }
1389*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,localX,&x));
1390*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localX));
1391e77315bcSPatrick Sanan 
1392*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops((41.0 + 8.0*POWFLOP)*xm*ym));
1393e77315bcSPatrick Sanan   PetscFunctionReturn(0);
1394e77315bcSPatrick Sanan }
1395e77315bcSPatrick Sanan 
1396e77315bcSPatrick Sanan /*TEST
1397e77315bcSPatrick Sanan 
1398e77315bcSPatrick Sanan    test:
1399e77315bcSPatrick Sanan       nsize: 4
1400e77315bcSPatrick 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
1401e77315bcSPatrick Sanan       requires: !single
1402e77315bcSPatrick Sanan 
1403e77315bcSPatrick Sanan TEST*/
1404