xref: /petsc/src/ts/tutorials/ex29.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Time-Dependent Allan-Cahn example in 2D with Varying Coefficients";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown /*
4*c4762a1bSJed Brown  This example is mainly here to show how to transfer coefficients between subdomains and levels in
5*c4762a1bSJed Brown  multigrid and domain decomposition.
6*c4762a1bSJed Brown  */
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown /*T
9*c4762a1bSJed Brown    Concepts: Alan-Cahn
10*c4762a1bSJed Brown    Concepts: DMDA with timestepping
11*c4762a1bSJed Brown    Concepts: Variable Coefficient
12*c4762a1bSJed Brown    Concepts: Nonlinear Domain Decomposition and Multigrid with Coefficients
13*c4762a1bSJed Brown    Processors: n
14*c4762a1bSJed Brown T*/
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown #include <petscdm.h>
17*c4762a1bSJed Brown #include <petscdmda.h>
18*c4762a1bSJed Brown #include <petscsnes.h>
19*c4762a1bSJed Brown #include <petscts.h>
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown typedef struct {
22*c4762a1bSJed Brown   PetscScalar epsilon;
23*c4762a1bSJed Brown   PetscScalar beta;
24*c4762a1bSJed Brown } Coeff;
25*c4762a1bSJed Brown 
26*c4762a1bSJed Brown typedef struct {
27*c4762a1bSJed Brown   PetscScalar u;
28*c4762a1bSJed Brown } Field;
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(DM da,void *ctx,Vec X);
31*c4762a1bSJed Brown extern PetscErrorCode FormDiffusionCoefficient(DM da,void *ctx,Vec X);
32*c4762a1bSJed Brown extern PetscErrorCode FormIFunctionLocal(DMDALocalInfo*,PetscReal,Field**,Field**,Field**,void*);
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown /* hooks */
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown static PetscErrorCode CoefficientCoarsenHook(DM dm, DM dmc,void *ctx)
37*c4762a1bSJed Brown {
38*c4762a1bSJed Brown   Vec            c,cc,ccl;
39*c4762a1bSJed Brown   PetscErrorCode ierr;
40*c4762a1bSJed Brown   Mat            J;
41*c4762a1bSJed Brown   Vec            vscale;
42*c4762a1bSJed Brown   DM             cdm,cdmc;
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown   PetscFunctionBegin;
45*c4762a1bSJed Brown   ierr = PetscObjectQuery((PetscObject)dm,"coefficientdm",(PetscObject*)&cdm);CHKERRQ(ierr);
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown   if (!cdm) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"The coefficient DM needs to be set up!");
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown   ierr = DMDACreateCompatibleDMDA(dmc,2,&cdmc);CHKERRQ(ierr);
50*c4762a1bSJed Brown   ierr = PetscObjectCompose((PetscObject)dmc,"coefficientdm",(PetscObject)cdmc);CHKERRQ(ierr);
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown   ierr = DMGetNamedGlobalVector(cdm,"coefficient",&c);CHKERRQ(ierr);
53*c4762a1bSJed Brown   ierr = DMGetNamedGlobalVector(cdmc,"coefficient",&cc);CHKERRQ(ierr);
54*c4762a1bSJed Brown   ierr = DMGetNamedLocalVector(cdmc,"coefficient",&ccl);CHKERRQ(ierr);
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown   ierr = DMCreateInterpolation(cdmc,cdm,&J,&vscale);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = MatRestrict(J,c,cc);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = VecPointwiseMult(cc,vscale,cc);CHKERRQ(ierr);
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = VecDestroy(&vscale);CHKERRQ(ierr);
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(cdmc,cc,INSERT_VALUES,ccl);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(cdmc,cc,INSERT_VALUES,ccl);CHKERRQ(ierr);
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown   ierr = DMRestoreNamedGlobalVector(cdm,"coefficient",&c);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = DMRestoreNamedGlobalVector(cdmc,"coefficient",&cc);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = DMRestoreNamedLocalVector(cdmc,"coefficient",&ccl);CHKERRQ(ierr);
69*c4762a1bSJed Brown 
70*c4762a1bSJed Brown   ierr = DMCoarsenHookAdd(dmc,CoefficientCoarsenHook,NULL,NULL);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = DMDestroy(&cdmc);CHKERRQ(ierr);
72*c4762a1bSJed Brown   PetscFunctionReturn(0);
73*c4762a1bSJed Brown }
74*c4762a1bSJed Brown 
75*c4762a1bSJed Brown /* This could restrict auxiliary information to the coarse level.
76*c4762a1bSJed Brown  */
77*c4762a1bSJed Brown static PetscErrorCode CoefficientSubDomainRestrictHook(DM dm,DM subdm,void *ctx)
78*c4762a1bSJed Brown {
79*c4762a1bSJed Brown   Vec            c,cc;
80*c4762a1bSJed Brown   DM             cdm,csubdm;
81*c4762a1bSJed Brown   PetscErrorCode ierr;
82*c4762a1bSJed Brown   VecScatter     *iscat,*oscat,*gscat;
83*c4762a1bSJed Brown 
84*c4762a1bSJed Brown   PetscFunctionBegin;
85*c4762a1bSJed Brown   ierr = PetscObjectQuery((PetscObject)dm,"coefficientdm",(PetscObject*)&cdm);CHKERRQ(ierr);
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown   if (!cdm) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"The coefficient DM needs to be set up!");
88*c4762a1bSJed Brown 
89*c4762a1bSJed Brown   ierr = DMDACreateCompatibleDMDA(subdm,2,&csubdm);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = PetscObjectCompose((PetscObject)subdm,"coefficientdm",(PetscObject)csubdm);CHKERRQ(ierr);
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   ierr = DMGetNamedGlobalVector(cdm,"coefficient",&c);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = DMGetNamedLocalVector(csubdm,"coefficient",&cc);CHKERRQ(ierr);
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   ierr = DMCreateDomainDecompositionScatters(cdm,1,&csubdm,&iscat,&oscat,&gscat);CHKERRQ(ierr);
96*c4762a1bSJed Brown 
97*c4762a1bSJed Brown   ierr = VecScatterBegin(*gscat,c,cc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = VecScatterEnd(*gscat,c,cc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   ierr = VecScatterDestroy(iscat);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = VecScatterDestroy(oscat);CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr = VecScatterDestroy(gscat);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr = PetscFree(iscat);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = PetscFree(oscat);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = PetscFree(gscat);CHKERRQ(ierr);
106*c4762a1bSJed Brown 
107*c4762a1bSJed Brown   ierr = DMRestoreNamedGlobalVector(cdm,"coefficient",&c);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = DMRestoreNamedLocalVector(csubdm,"coefficient",&cc);CHKERRQ(ierr);
109*c4762a1bSJed Brown 
110*c4762a1bSJed Brown   ierr = DMDestroy(&csubdm);CHKERRQ(ierr);
111*c4762a1bSJed Brown   PetscFunctionReturn(0);
112*c4762a1bSJed Brown }
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown int main(int argc,char **argv)
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown {
117*c4762a1bSJed Brown   TS             ts;
118*c4762a1bSJed Brown   Vec            x,c,clocal;
119*c4762a1bSJed Brown   PetscErrorCode ierr;
120*c4762a1bSJed Brown   DM             da,cda;
121*c4762a1bSJed Brown 
122*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
123*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
126*c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr);
127*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
128*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
129*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr);
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown   ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr);
132*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr);
133*c4762a1bSJed Brown 
134*c4762a1bSJed Brown   ierr = TSSetDM(ts, da);CHKERRQ(ierr);
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown   ierr = FormInitialGuess(da,NULL,x);CHKERRQ(ierr);
137*c4762a1bSJed Brown   ierr = DMDATSSetIFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*))FormIFunctionLocal,NULL);CHKERRQ(ierr);
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown   /* set up the coefficient */
140*c4762a1bSJed Brown 
141*c4762a1bSJed Brown   ierr = DMDACreateCompatibleDMDA(da,2,&cda);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = PetscObjectCompose((PetscObject)da,"coefficientdm",(PetscObject)cda);CHKERRQ(ierr);
143*c4762a1bSJed Brown 
144*c4762a1bSJed Brown   ierr = DMGetNamedGlobalVector(cda,"coefficient",&c);CHKERRQ(ierr);
145*c4762a1bSJed Brown   ierr = DMGetNamedLocalVector(cda,"coefficient",&clocal);CHKERRQ(ierr);
146*c4762a1bSJed Brown 
147*c4762a1bSJed Brown   ierr = FormDiffusionCoefficient(cda,NULL,c);CHKERRQ(ierr);
148*c4762a1bSJed Brown 
149*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(cda,c,INSERT_VALUES,clocal);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(cda,c,INSERT_VALUES,clocal);CHKERRQ(ierr);
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown   ierr = DMRestoreNamedLocalVector(cda,"coefficient",&clocal);CHKERRQ(ierr);
153*c4762a1bSJed Brown   ierr = DMRestoreNamedGlobalVector(cda,"coefficient",&c);CHKERRQ(ierr);
154*c4762a1bSJed Brown 
155*c4762a1bSJed Brown   ierr = DMCoarsenHookAdd(da,CoefficientCoarsenHook,NULL,NULL);CHKERRQ(ierr);
156*c4762a1bSJed Brown   ierr = DMSubDomainHookAdd(da,CoefficientSubDomainRestrictHook,NULL,NULL);CHKERRQ(ierr);
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,10000);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,10000.0);CHKERRQ(ierr);
160*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
161*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,0.05);CHKERRQ(ierr);
162*c4762a1bSJed Brown   ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
163*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown   ierr = TSSolve(ts,x);CHKERRQ(ierr);
166*c4762a1bSJed Brown 
167*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
169*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
170*c4762a1bSJed Brown   ierr = DMDestroy(&cda);CHKERRQ(ierr);
171*c4762a1bSJed Brown 
172*c4762a1bSJed Brown   ierr = PetscFinalize();
173*c4762a1bSJed Brown   return ierr;
174*c4762a1bSJed Brown }
175*c4762a1bSJed Brown 
176*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
177*c4762a1bSJed Brown 
178*c4762a1bSJed Brown PetscErrorCode FormInitialGuess(DM da,void *ctx,Vec X)
179*c4762a1bSJed Brown {
180*c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
181*c4762a1bSJed Brown   PetscErrorCode ierr;
182*c4762a1bSJed Brown   Field          **x;
183*c4762a1bSJed Brown   PetscReal      x0,x1;
184*c4762a1bSJed Brown 
185*c4762a1bSJed Brown   PetscFunctionBeginUser;
186*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr);
189*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
190*c4762a1bSJed Brown 
191*c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
192*c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
193*c4762a1bSJed Brown       x0        = 10.0*(i - 0.5*(Mx-1)) / (Mx-1);
194*c4762a1bSJed Brown       x1        = 10.0*(j - 0.5*(Mx-1)) / (My-1);
195*c4762a1bSJed Brown       x[j][i].u = PetscCosReal(2.0*PetscSqrtReal(x1*x1 + x0*x0));
196*c4762a1bSJed Brown     }
197*c4762a1bSJed Brown   }
198*c4762a1bSJed Brown 
199*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr);
200*c4762a1bSJed Brown   PetscFunctionReturn(0);
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown }
203*c4762a1bSJed Brown 
204*c4762a1bSJed Brown PetscErrorCode FormDiffusionCoefficient(DM da,void *ctx,Vec X)
205*c4762a1bSJed Brown {
206*c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
207*c4762a1bSJed Brown   PetscErrorCode ierr;
208*c4762a1bSJed Brown   Coeff          **x;
209*c4762a1bSJed Brown   PetscReal      x1,x0;
210*c4762a1bSJed Brown 
211*c4762a1bSJed Brown   PetscFunctionBeginUser;
212*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
213*c4762a1bSJed Brown 
214*c4762a1bSJed Brown   /*
215*c4762a1bSJed Brown   ierr = VecSetRandom(X,NULL);
216*c4762a1bSJed Brown   ierr = VecMin(X,NULL,&min);CHKERRQ(ierr);
217*c4762a1bSJed Brown    */
218*c4762a1bSJed Brown 
219*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr);
220*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
221*c4762a1bSJed Brown 
222*c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
223*c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
224*c4762a1bSJed Brown       x0 = 10.0*(i - 0.5*(Mx-1)) / (Mx-1);
225*c4762a1bSJed Brown       x1 = 10.0*(j - 0.5*(My-1)) / (My-1);
226*c4762a1bSJed Brown 
227*c4762a1bSJed Brown       x[j][i].epsilon = 0.0;
228*c4762a1bSJed Brown       x[j][i].beta    = 0.05+0.05*PetscSqrtReal(x0*x0+x1*x1);
229*c4762a1bSJed Brown     }
230*c4762a1bSJed Brown   }
231*c4762a1bSJed Brown 
232*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr);
233*c4762a1bSJed Brown   PetscFunctionReturn(0);
234*c4762a1bSJed Brown 
235*c4762a1bSJed Brown }
236*c4762a1bSJed Brown 
237*c4762a1bSJed Brown PetscErrorCode FormIFunctionLocal(DMDALocalInfo *info,PetscReal ptime,Field **x,Field **xt,Field **f,void *ctx)
238*c4762a1bSJed Brown {
239*c4762a1bSJed Brown   PetscErrorCode ierr;
240*c4762a1bSJed Brown   PetscInt       i,j;
241*c4762a1bSJed Brown   PetscReal      hx,hy,dhx,dhy,hxdhy,hydhx,scale;
242*c4762a1bSJed Brown   PetscScalar    u,uxx,uyy;
243*c4762a1bSJed Brown   PetscScalar    ux,uy,bx,by;
244*c4762a1bSJed Brown   Vec            C;
245*c4762a1bSJed Brown   Coeff          **c;
246*c4762a1bSJed Brown   DM             cdm;
247*c4762a1bSJed Brown 
248*c4762a1bSJed Brown   PetscFunctionBeginUser;
249*c4762a1bSJed Brown   ierr = PetscObjectQuery((PetscObject)info->da,"coefficientdm",(PetscObject*)&cdm);CHKERRQ(ierr);
250*c4762a1bSJed Brown   ierr = DMGetNamedLocalVector(cdm,"coefficient",&C);CHKERRQ(ierr);
251*c4762a1bSJed Brown   ierr = DMDAVecGetArray(cdm,C,&c);CHKERRQ(ierr);
252*c4762a1bSJed Brown 
253*c4762a1bSJed Brown   hx = 10.0/((PetscReal)(info->mx-1));
254*c4762a1bSJed Brown   hy = 10.0/((PetscReal)(info->my-1));
255*c4762a1bSJed Brown 
256*c4762a1bSJed Brown   dhx = 1. / hx;
257*c4762a1bSJed Brown   dhy = 1. / hy;
258*c4762a1bSJed Brown 
259*c4762a1bSJed Brown   hxdhy =  hx/hy;
260*c4762a1bSJed Brown   hydhx =  hy/hx;
261*c4762a1bSJed Brown   scale =   hx*hy;
262*c4762a1bSJed Brown 
263*c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
264*c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
265*c4762a1bSJed Brown       f[j][i].u = xt[j][i].u*scale;
266*c4762a1bSJed Brown 
267*c4762a1bSJed Brown       u = x[j][i].u;
268*c4762a1bSJed Brown 
269*c4762a1bSJed Brown       f[j][i].u += scale*(u*u - 1.)*u;
270*c4762a1bSJed Brown 
271*c4762a1bSJed Brown       if (i == 0)               f[j][i].u += (x[j][i].u - x[j][i+1].u)*dhx;
272*c4762a1bSJed Brown       else if (i == info->mx-1) f[j][i].u += (x[j][i].u - x[j][i-1].u)*dhx;
273*c4762a1bSJed Brown       else if (j == 0)          f[j][i].u += (x[j][i].u - x[j+1][i].u)*dhy;
274*c4762a1bSJed Brown       else if (j == info->my-1) f[j][i].u += (x[j][i].u - x[j-1][i].u)*dhy;
275*c4762a1bSJed Brown       else {
276*c4762a1bSJed Brown         uyy     = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
277*c4762a1bSJed Brown         uxx     = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
278*c4762a1bSJed Brown 
279*c4762a1bSJed Brown         bx      = 0.5*(c[j][i+1].beta - c[j][i-1].beta)*dhx;
280*c4762a1bSJed Brown         by      = 0.5*(c[j+1][i].beta - c[j-1][i].beta)*dhy;
281*c4762a1bSJed Brown 
282*c4762a1bSJed Brown         ux      = 0.5*(x[j][i+1].u - x[j][i-1].u)*dhx;
283*c4762a1bSJed Brown         uy      = 0.5*(x[j+1][i].u - x[j-1][i].u)*dhy;
284*c4762a1bSJed Brown 
285*c4762a1bSJed Brown         f[j][i].u += c[j][i].beta*(uxx + uyy) + scale*(bx*ux + by*uy);
286*c4762a1bSJed Brown        }
287*c4762a1bSJed Brown     }
288*c4762a1bSJed Brown   }
289*c4762a1bSJed Brown   ierr = PetscLogFlops(11.*info->ym*info->xm);CHKERRQ(ierr);
290*c4762a1bSJed Brown 
291*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(cdm,C,&c);CHKERRQ(ierr);
292*c4762a1bSJed Brown   ierr = DMRestoreNamedLocalVector(cdm,"coefficient",&C);CHKERRQ(ierr);
293*c4762a1bSJed Brown   PetscFunctionReturn(0);
294*c4762a1bSJed Brown }
295*c4762a1bSJed Brown 
296*c4762a1bSJed Brown /*TEST
297*c4762a1bSJed Brown 
298*c4762a1bSJed Brown     test:
299*c4762a1bSJed Brown       args: -da_refine 4 -ts_max_steps 10 -ts_rtol 1e-3 -ts_atol 1e-3 -ts_type arkimex -ts_monitor -snes_monitor -snes_type ngmres  -npc_snes_type nasm -npc_snes_nasm_type restrict -da_overlap 4
300*c4762a1bSJed Brown       nsize: 16
301*c4762a1bSJed Brown       requires: !single
302*c4762a1bSJed Brown       output_file: output/ex29.out
303*c4762a1bSJed Brown 
304*c4762a1bSJed Brown TEST*/
305