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