Actual source code: ex18.c
1: /* $Id: ex18.c,v 1.23 2001/08/07 21:31:17 bsmith Exp $ */
4: static char help[] ="Nonlinear Radiative Transport PDE with multigrid in 2d.n
5: Uses 2-dimensional distributed arrays.n
6: A 2-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 = 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 5-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 David Keyes
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,3,&user,&dmmg);
88: /*
89: Set the DA (grid structure) for the grids.
90: */
91: DACreate2d(PETSC_COMM_WORLD,DA_NONPERIODIC,DA_STENCIL_STAR,5,5,PETSC_DECIDE,PETSC_DECIDE,1,1,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,ierr,xs,ys,xm,ym;
130: PetscReal tleft = user->tleft;
131: PetscScalar **x;
135: /* Get ghost points */
136: DAGetCorners((DA)dmmg->dm,&xs,&ys,0,&xm,&ym,0);
137: DAVecGetArray((DA)dmmg->dm,X,(void**)&x);
139: /* Compute initial guess */
140: for (j=ys; j<ys+ym; j++) {
141: for (i=xs; i<xs+xm; i++) {
142: x[j][i] = tleft;
143: }
144: }
145: DAVecRestoreArray((DA)dmmg->dm,X,(void**)&x);
146: return(0);
147: }
148: /* -------------------- Evaluate Function F(x) --------------------- */
149: int FormFunction(SNES snes,Vec X,Vec F,void* ptr)
150: {
151: DMMG dmmg = (DMMG)ptr;
152: AppCtx *user = (AppCtx*)dmmg->user;
153: int ierr,i,j,mx,my,xs,ys,xm,ym;
154: PetscScalar zero = 0.0,one = 1.0;
155: PetscScalar hx,hy,hxdhy,hydhx;
156: 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;
157: PetscScalar tleft,tright,beta;
158: PetscScalar **x,**f;
159: Vec localX;
162: DAGetLocalVector((DA)dmmg->dm,&localX);
163: DAGetInfo((DA)dmmg->dm,PETSC_NULL,&mx,&my,0,0,0,0,0,0,0,0);
164: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
165: hxdhy = hx/hy; hydhx = hy/hx;
166: tleft = user->tleft; tright = user->tright;
167: beta = user->beta;
168:
169: /* Get ghost points */
170: DAGlobalToLocalBegin((DA)dmmg->dm,X,INSERT_VALUES,localX);
171: DAGlobalToLocalEnd((DA)dmmg->dm,X,INSERT_VALUES,localX);
172: DAGetCorners((DA)dmmg->dm,&xs,&ys,0,&xm,&ym,0);
173: DAVecGetArray((DA)dmmg->dm,localX,(void**)&x);
174: DAVecGetArray((DA)dmmg->dm,F,(void**)&f);
176: /* Evaluate function */
177: for (j=ys; j<ys+ym; j++) {
178: for (i=xs; i<xs+xm; i++) {
179: t0 = x[j][i];
181: if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
183: /* general interior volume */
185: tw = x[j][i-1];
186: aw = 0.5*(t0 + tw);
187: dw = PetscPowScalar(aw,beta);
188: fw = dw*(t0 - tw);
190: te = x[j][i+1];
191: ae = 0.5*(t0 + te);
192: de = PetscPowScalar(ae,beta);
193: fe = de*(te - t0);
195: ts = x[j-1][i];
196: as = 0.5*(t0 + ts);
197: ds = PetscPowScalar(as,beta);
198: fs = ds*(t0 - ts);
199:
200: tn = x[j+1][i];
201: an = 0.5*(t0 + tn);
202: dn = PetscPowScalar(an,beta);
203: fn = dn*(tn - t0);
205: } else if (i == 0) {
207: /* left-hand boundary */
208: tw = tleft;
209: aw = 0.5*(t0 + tw);
210: dw = PetscPowScalar(aw,beta);
211: fw = dw*(t0 - tw);
213: te = x[j][i+1];
214: ae = 0.5*(t0 + te);
215: de = PetscPowScalar(ae,beta);
216: fe = de*(te - t0);
218: if (j > 0) {
219: ts = x[j-1][i];
220: as = 0.5*(t0 + ts);
221: ds = PetscPowScalar(as,beta);
222: fs = ds*(t0 - ts);
223: } else {
224: fs = zero;
225: }
227: if (j < my-1) {
228: tn = x[j+1][i];
229: an = 0.5*(t0 + tn);
230: dn = PetscPowScalar(an,beta);
231: fn = dn*(tn - t0);
232: } else {
233: fn = zero;
234: }
236: } else if (i == mx-1) {
238: /* right-hand boundary */
239: tw = x[j][i-1];
240: aw = 0.5*(t0 + tw);
241: dw = PetscPowScalar(aw,beta);
242: fw = dw*(t0 - tw);
243:
244: te = tright;
245: ae = 0.5*(t0 + te);
246: de = PetscPowScalar(ae,beta);
247: fe = de*(te - t0);
248:
249: if (j > 0) {
250: ts = x[j-1][i];
251: as = 0.5*(t0 + ts);
252: ds = PetscPowScalar(as,beta);
253: fs = ds*(t0 - ts);
254: } else {
255: fs = zero;
256: }
257:
258: if (j < my-1) {
259: tn = x[j+1][i];
260: an = 0.5*(t0 + tn);
261: dn = PetscPowScalar(an,beta);
262: fn = dn*(tn - t0);
263: } else {
264: fn = zero;
265: }
267: } else if (j == 0) {
269: /* bottom boundary,and i <> 0 or mx-1 */
270: tw = x[j][i-1];
271: aw = 0.5*(t0 + tw);
272: dw = PetscPowScalar(aw,beta);
273: fw = dw*(t0 - tw);
275: te = x[j][i+1];
276: ae = 0.5*(t0 + te);
277: de = PetscPowScalar(ae,beta);
278: fe = de*(te - t0);
280: fs = zero;
282: tn = x[j+1][i];
283: an = 0.5*(t0 + tn);
284: dn = PetscPowScalar(an,beta);
285: fn = dn*(tn - t0);
287: } else if (j == my-1) {
289: /* top boundary,and i <> 0 or mx-1 */
290: tw = x[j][i-1];
291: aw = 0.5*(t0 + tw);
292: dw = PetscPowScalar(aw,beta);
293: fw = dw*(t0 - tw);
295: te = x[j][i+1];
296: ae = 0.5*(t0 + te);
297: de = PetscPowScalar(ae,beta);
298: fe = de*(te - t0);
300: ts = x[j-1][i];
301: as = 0.5*(t0 + ts);
302: ds = PetscPowScalar(as,beta);
303: fs = ds*(t0 - ts);
305: fn = zero;
307: }
309: f[j][i] = - hydhx*(fe-fw) - hxdhy*(fn-fs);
311: }
312: }
313: DAVecRestoreArray((DA)dmmg->dm,localX,(void**)&x);
314: DAVecRestoreArray((DA)dmmg->dm,F,(void**)&f);
315: DARestoreLocalVector((DA)dmmg->dm,&localX);
316: PetscLogFlops((22 + 4*POWFLOP)*ym*xm);
317: return(0);
318: }
319: /* -------------------- Evaluate Jacobian F(x) --------------------- */
320: int FormJacobian(SNES snes,Vec X,Mat *J,Mat *B,MatStructure *flg,void *ptr)
321: {
322: DMMG dmmg = (DMMG)ptr;
323: AppCtx *user = (AppCtx*)dmmg->user;
324: Mat jac = *J;
325: int ierr,i,j,mx,my,xs,ys,xm,ym;
326: PetscScalar one = 1.0,hx,hy,hxdhy,hydhx,t0,tn,ts,te,tw;
327: PetscScalar dn,ds,de,dw,an,as,ae,aw,bn,bs,be,bw,gn,gs,ge,gw;
328: PetscScalar tleft,tright,beta,bm1,coef;
329: PetscScalar v[5],**x;
330: Vec localX;
331: MatStencil col[5],row;
334: DAGetLocalVector((DA)dmmg->dm,&localX);
335: *flg = SAME_NONZERO_PATTERN;
336: DAGetInfo((DA)dmmg->dm,PETSC_NULL,&mx,&my,0,0,0,0,0,0,0,0);
337: hx = one/(PetscReal)(mx-1); hy = one/(PetscReal)(my-1);
338: hxdhy = hx/hy; hydhx = hy/hx;
339: tleft = user->tleft; tright = user->tright;
340: beta = user->beta; bm1 = user->bm1; coef = user->coef;
342: /* Get ghost points */
343: DAGlobalToLocalBegin((DA)dmmg->dm,X,INSERT_VALUES,localX);
344: DAGlobalToLocalEnd((DA)dmmg->dm,X,INSERT_VALUES,localX);
345: DAGetCorners((DA)dmmg->dm,&xs,&ys,0,&xm,&ym,0);
346: DAVecGetArray((DA)dmmg->dm,localX,(void**)&x);
348: /* Evaluate Jacobian of function */
349: for (j=ys; j<ys+ym; j++) {
350: for (i=xs; i<xs+xm; i++) {
351: t0 = x[j][i];
353: if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
355: /* general interior volume */
357: tw = x[j][i-1];
358: aw = 0.5*(t0 + tw);
359: bw = PetscPowScalar(aw,bm1);
360: /* dw = bw * aw */
361: dw = PetscPowScalar(aw,beta);
362: gw = coef*bw*(t0 - tw);
364: te = x[j][i+1];
365: ae = 0.5*(t0 + te);
366: be = PetscPowScalar(ae,bm1);
367: /* de = be * ae; */
368: de = PetscPowScalar(ae,beta);
369: ge = coef*be*(te - t0);
371: ts = x[j-1][i];
372: as = 0.5*(t0 + ts);
373: bs = PetscPowScalar(as,bm1);
374: /* ds = bs * as; */
375: ds = PetscPowScalar(as,beta);
376: gs = coef*bs*(t0 - ts);
377:
378: tn = x[j+1][i];
379: an = 0.5*(t0 + tn);
380: bn = PetscPowScalar(an,bm1);
381: /* dn = bn * an; */
382: dn = PetscPowScalar(an,beta);
383: gn = coef*bn*(tn - t0);
385: v[0] = - hxdhy*(ds - gs); col[0].j = j-1; col[0].i = i;
386: v[1] = - hydhx*(dw - gw); col[1].j = j; col[1].i = i-1;
387: v[2] = hxdhy*(ds + dn + gs - gn) + hydhx*(dw + de + gw - ge); col[2].j = row.j = j; col[2].i = row.i = i;
388: v[3] = - hydhx*(de + ge); col[3].j = j; col[3].i = i+1;
389: v[4] = - hxdhy*(dn + gn); col[4].j = j+1; col[4].i = i;
390: MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
392: } else if (i == 0) {
394: /* left-hand boundary */
395: tw = tleft;
396: aw = 0.5*(t0 + tw);
397: bw = PetscPowScalar(aw,bm1);
398: /* dw = bw * aw */
399: dw = PetscPowScalar(aw,beta);
400: gw = coef*bw*(t0 - tw);
401:
402: te = x[j][i + 1];
403: ae = 0.5*(t0 + te);
404: be = PetscPowScalar(ae,bm1);
405: /* de = be * ae; */
406: de = PetscPowScalar(ae,beta);
407: ge = coef*be*(te - t0);
408:
409: /* left-hand bottom boundary */
410: if (j == 0) {
412: tn = x[j+1][i];
413: an = 0.5*(t0 + tn);
414: bn = PetscPowScalar(an,bm1);
415: /* dn = bn * an; */
416: dn = PetscPowScalar(an,beta);
417: gn = coef*bn*(tn - t0);
418:
419: v[0] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[0].j = row.j = j; col[0].i = row.i = i;
420: v[1] = - hydhx*(de + ge); col[1].j = j; col[1].i = i+1;
421: v[2] = - hxdhy*(dn + gn); col[2].j = j+1; col[2].i = i;
422: MatSetValuesStencil(jac,1,&row,3,col,v,INSERT_VALUES);
423:
424: /* left-hand interior boundary */
425: } else if (j < my-1) {
427: ts = x[j-1][i];
428: as = 0.5*(t0 + ts);
429: bs = PetscPowScalar(as,bm1);
430: /* ds = bs * as; */
431: ds = PetscPowScalar(as,beta);
432: gs = coef*bs*(t0 - ts);
433:
434: tn = x[j+1][i];
435: an = 0.5*(t0 + tn);
436: bn = PetscPowScalar(an,bm1);
437: /* dn = bn * an; */
438: dn = PetscPowScalar(an,beta);
439: gn = coef*bn*(tn - t0);
440:
441: v[0] = - hxdhy*(ds - gs); col[0].j = j-1; col[0].i = i;
442: v[1] = hxdhy*(ds + dn + gs - gn) + hydhx*(dw + de + gw - ge); col[1].j = row.j = j; col[1].i = row.i = i;
443: v[2] = - hydhx*(de + ge); col[2].j = j; col[2].i = i+1;
444: v[3] = - hxdhy*(dn + gn); col[3].j = j+1; col[3].i = i;
445: MatSetValuesStencil(jac,1,&row,4,col,v,INSERT_VALUES);
446: /* left-hand top boundary */
447: } else {
449: ts = x[j-1][i];
450: as = 0.5*(t0 + ts);
451: bs = PetscPowScalar(as,bm1);
452: /* ds = bs * as; */
453: ds = PetscPowScalar(as,beta);
454: gs = coef*bs*(t0 - ts);
455:
456: v[0] = - hxdhy*(ds - gs); col[0].j = j-1; col[0].i = i;
457: v[1] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge); col[1].j = row.j = j; col[1].i = row.i = i;
458: v[2] = - hydhx*(de + ge); col[2].j = j; col[2].i = i+1;
459: MatSetValuesStencil(jac,1,&row,3,col,v,INSERT_VALUES);
460: }
462: } else if (i == mx-1) {
463:
464: /* right-hand boundary */
465: tw = x[j][i-1];
466: aw = 0.5*(t0 + tw);
467: bw = PetscPowScalar(aw,bm1);
468: /* dw = bw * aw */
469: dw = PetscPowScalar(aw,beta);
470: gw = coef*bw*(t0 - tw);
471:
472: te = tright;
473: ae = 0.5*(t0 + te);
474: be = PetscPowScalar(ae,bm1);
475: /* de = be * ae; */
476: de = PetscPowScalar(ae,beta);
477: ge = coef*be*(te - t0);
478:
479: /* right-hand bottom boundary */
480: if (j == 0) {
482: tn = x[j+1][i];
483: an = 0.5*(t0 + tn);
484: bn = PetscPowScalar(an,bm1);
485: /* dn = bn * an; */
486: dn = PetscPowScalar(an,beta);
487: gn = coef*bn*(tn - t0);
488:
489: v[0] = - hydhx*(dw - gw); col[0].j = j; col[0].i = i-1;
490: v[1] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[1].j = row.j = j; col[1].i = row.i = i;
491: v[2] = - hxdhy*(dn + gn); col[2].j = j+1; col[2].i = i;
492: MatSetValuesStencil(jac,1,&row,3,col,v,INSERT_VALUES);
493:
494: /* right-hand interior boundary */
495: } else if (j < my-1) {
497: ts = x[j-1][i];
498: as = 0.5*(t0 + ts);
499: bs = PetscPowScalar(as,bm1);
500: /* ds = bs * as; */
501: ds = PetscPowScalar(as,beta);
502: gs = coef*bs*(t0 - ts);
503:
504: tn = x[j+1][i];
505: an = 0.5*(t0 + tn);
506: bn = PetscPowScalar(an,bm1);
507: /* dn = bn * an; */
508: dn = PetscPowScalar(an,beta);
509: gn = coef*bn*(tn - t0);
510:
511: v[0] = - hxdhy*(ds - gs); col[0].j = j-1; col[0].i = i;
512: v[1] = - hydhx*(dw - gw); col[1].j = j; col[1].i = i-1;
513: v[2] = hxdhy*(ds + dn + gs - gn) + hydhx*(dw + de + gw - ge); col[2].j = row.j = j; col[2].i = row.i = i;
514: v[3] = - hxdhy*(dn + gn); col[3].j = j+1; col[3].i = i;
515: MatSetValuesStencil(jac,1,&row,4,col,v,INSERT_VALUES);
516: /* right-hand top boundary */
517: } else {
519: ts = x[j-1][i];
520: as = 0.5*(t0 + ts);
521: bs = PetscPowScalar(as,bm1);
522: /* ds = bs * as; */
523: ds = PetscPowScalar(as,beta);
524: gs = coef*bs*(t0 - ts);
525:
526: v[0] = - hxdhy*(ds - gs); col[0].j = j-1; col[0].i = i;
527: v[1] = - hydhx*(dw - gw); col[1].j = j; col[1].i = i-1;
528: v[2] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge); col[2].j = row.j = j; col[2].i = row.i = i;
529: MatSetValuesStencil(jac,1,&row,3,col,v,INSERT_VALUES);
530: }
532: /* bottom boundary,and i <> 0 or mx-1 */
533: } else if (j == 0) {
535: tw = x[j][i-1];
536: aw = 0.5*(t0 + tw);
537: bw = PetscPowScalar(aw,bm1);
538: /* dw = bw * aw */
539: dw = PetscPowScalar(aw,beta);
540: gw = coef*bw*(t0 - tw);
542: te = x[j][i+1];
543: ae = 0.5*(t0 + te);
544: be = PetscPowScalar(ae,bm1);
545: /* de = be * ae; */
546: de = PetscPowScalar(ae,beta);
547: ge = coef*be*(te - t0);
549: tn = x[j+1][i];
550: an = 0.5*(t0 + tn);
551: bn = PetscPowScalar(an,bm1);
552: /* dn = bn * an; */
553: dn = PetscPowScalar(an,beta);
554: gn = coef*bn*(tn - t0);
555:
556: v[0] = - hydhx*(dw - gw); col[0].j = j; col[0].i = i-1;
557: v[1] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[1].j = row.j = j; col[1].i = row.i = i;
558: v[2] = - hydhx*(de + ge); col[2].j = j; col[2].i = i+1;
559: v[3] = - hxdhy*(dn + gn); col[3].j = j+1; col[3].i = i;
560: MatSetValuesStencil(jac,1,&row,4,col,v,INSERT_VALUES);
561:
562: /* top boundary,and i <> 0 or mx-1 */
563: } else if (j == my-1) {
564:
565: tw = x[j][i-1];
566: aw = 0.5*(t0 + tw);
567: bw = PetscPowScalar(aw,bm1);
568: /* dw = bw * aw */
569: dw = PetscPowScalar(aw,beta);
570: gw = coef*bw*(t0 - tw);
572: te = x[j][i+1];
573: ae = 0.5*(t0 + te);
574: be = PetscPowScalar(ae,bm1);
575: /* de = be * ae; */
576: de = PetscPowScalar(ae,beta);
577: ge = coef*be*(te - t0);
579: ts = x[j-1][i];
580: as = 0.5*(t0 + ts);
581: bs = PetscPowScalar(as,bm1);
582: /* ds = bs * as; */
583: ds = PetscPowScalar(as,beta);
584: gs = coef*bs*(t0 - ts);
586: v[0] = - hxdhy*(ds - gs); col[0].j = j-1; col[0].i = i;
587: v[1] = - hydhx*(dw - gw); col[1].j = j; col[1].i = i-1;
588: v[2] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge); col[2].j = row.j = j; col[2].i = row.i = i;
589: v[3] = - hydhx*(de + ge); col[3].j = j; col[3].i = i+1;
590: MatSetValuesStencil(jac,1,&row,4,col,v,INSERT_VALUES);
591:
592: }
593: }
594: }
595: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
596: DAVecRestoreArray((DA)dmmg->dm,localX,(void**)&x);
597: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
598: DARestoreLocalVector((DA)dmmg->dm,&localX);
600: PetscLogFlops((41 + 8*POWFLOP)*xm*ym);
601: return(0);
602: }