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