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