xref: /petsc/src/snes/tutorials/ex18.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] ="Nonlinear Radiative Transport PDE with multigrid in 2d.\n\
3c4762a1bSJed Brown Uses 2-dimensional distributed arrays.\n\
4c4762a1bSJed Brown A 2-dim simplified Radiative Transport test problem is used, with analytic Jacobian. \n\
5c4762a1bSJed Brown \n\
6c4762a1bSJed Brown   Solves the linear systems via multilevel methods \n\
7c4762a1bSJed Brown \n\
8c4762a1bSJed Brown The command line\n\
9c4762a1bSJed Brown options are:\n\
10c4762a1bSJed Brown   -tleft <tl>, where <tl> indicates the left Diriclet BC \n\
11c4762a1bSJed Brown   -tright <tr>, where <tr> indicates the right Diriclet BC \n\
12c4762a1bSJed Brown   -beta <beta>, where <beta> indicates the exponent in T \n\n";
13c4762a1bSJed Brown 
14c4762a1bSJed Brown /*T
15c4762a1bSJed Brown    Concepts: SNES^solving a system of nonlinear equations
16c4762a1bSJed Brown    Concepts: DMDA^using distributed arrays
17c4762a1bSJed Brown    Concepts: multigrid;
18c4762a1bSJed Brown    Processors: n
19c4762a1bSJed Brown T*/
20c4762a1bSJed Brown 
21c4762a1bSJed Brown /*
22c4762a1bSJed Brown 
23c4762a1bSJed Brown     This example models the partial differential equation
24c4762a1bSJed Brown 
25c4762a1bSJed Brown          - Div(alpha* T^beta (GRAD T)) = 0.
26c4762a1bSJed Brown 
27c4762a1bSJed Brown     where beta = 2.5 and alpha = 1.0
28c4762a1bSJed Brown 
29c4762a1bSJed Brown     BC: T_left = 1.0, T_right = 0.1, dT/dn_top = dTdn_bottom = 0.
30c4762a1bSJed Brown 
31c4762a1bSJed Brown     in the unit square, which is uniformly discretized in each of x and
32c4762a1bSJed Brown     y in this simple encoding.  The degrees of freedom are cell centered.
33c4762a1bSJed Brown 
34c4762a1bSJed Brown     A finite volume approximation with the usual 5-point stencil
35c4762a1bSJed Brown     is used to discretize the boundary value problem to obtain a
36c4762a1bSJed Brown     nonlinear system of equations.
37c4762a1bSJed Brown 
38c4762a1bSJed Brown     This code was contributed by David Keyes
39c4762a1bSJed Brown 
40c4762a1bSJed Brown */
41c4762a1bSJed Brown 
42c4762a1bSJed Brown #include <petscsnes.h>
43c4762a1bSJed Brown #include <petscdm.h>
44c4762a1bSJed Brown #include <petscdmda.h>
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /* User-defined application context */
47c4762a1bSJed Brown 
48c4762a1bSJed Brown typedef struct {
49c4762a1bSJed Brown   PetscReal tleft,tright;    /* Dirichlet boundary conditions */
50c4762a1bSJed Brown   PetscReal beta,bm1,coef;   /* nonlinear diffusivity parameterizations */
51c4762a1bSJed Brown } AppCtx;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown #define POWFLOP 5 /* assume a pow() takes five flops */
54c4762a1bSJed Brown 
55c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(SNES,Vec,void*);
56c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
57c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
58c4762a1bSJed Brown 
59c4762a1bSJed Brown int main(int argc,char **argv)
60c4762a1bSJed Brown {
61c4762a1bSJed Brown   SNES           snes;
62c4762a1bSJed Brown   AppCtx         user;
63c4762a1bSJed Brown   PetscInt       its,lits;
64c4762a1bSJed Brown   PetscReal      litspit;
65c4762a1bSJed Brown   DM             da;
66c4762a1bSJed Brown 
67*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* set problem parameters */
70c4762a1bSJed Brown   user.tleft  = 1.0;
71c4762a1bSJed Brown   user.tright = 0.1;
72c4762a1bSJed Brown   user.beta   = 2.5;
735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tleft",&user.tleft,NULL));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tright",&user.tright,NULL));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-beta",&user.beta,NULL));
76c4762a1bSJed Brown   user.bm1    = user.beta - 1.0;
77c4762a1bSJed Brown   user.coef   = user.beta/2.0;
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /*
80c4762a1bSJed Brown       Create the multilevel DM data structure
81c4762a1bSJed Brown   */
825f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /*
85c4762a1bSJed Brown       Set the DMDA (grid structure) for the grids.
86c4762a1bSJed Brown   */
875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,5,5,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(da,&user));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes,(DM)da));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /*
94c4762a1bSJed Brown      Create the nonlinear solver, and tell it the functions to use
95c4762a1bSJed Brown   */
965f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,NULL,FormFunction,&user));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,NULL,NULL,FormJacobian,&user));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetComputeInitialGuess(snes,FormInitialGuess,NULL));
100c4762a1bSJed Brown 
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,NULL,NULL));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetLinearSolveIterations(snes,&lits));
104c4762a1bSJed Brown   litspit = ((PetscReal)lits)/((PetscReal)its);
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %D\n",its));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of Linear iterations = %D\n",lits));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Average Linear its / SNES = %e\n",(double)litspit));
108c4762a1bSJed Brown 
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
111*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
112*b122ec5aSJacob Faibussowitsch   return 0;
113c4762a1bSJed Brown }
114c4762a1bSJed Brown /* --------------------  Form initial approximation ----------------- */
115c4762a1bSJed Brown PetscErrorCode FormInitialGuess(SNES snes,Vec X,void *ctx)
116c4762a1bSJed Brown {
117c4762a1bSJed Brown   AppCtx         *user;
118c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym;
119c4762a1bSJed Brown   PetscReal      tleft;
120c4762a1bSJed Brown   PetscScalar    **x;
121c4762a1bSJed Brown   DM             da;
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   PetscFunctionBeginUser;
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(da,&user));
126c4762a1bSJed Brown   tleft = user->tleft;
127c4762a1bSJed Brown   /* Get ghost points */
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /* Compute initial guess */
132c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
133c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
134c4762a1bSJed Brown       x[j][i] = tleft;
135c4762a1bSJed Brown     }
136c4762a1bSJed Brown   }
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
138c4762a1bSJed Brown   PetscFunctionReturn(0);
139c4762a1bSJed Brown }
140c4762a1bSJed Brown /* --------------------  Evaluate Function F(x) --------------------- */
141c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
142c4762a1bSJed Brown {
143c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
144c4762a1bSJed Brown   PetscInt       i,j,mx,my,xs,ys,xm,ym;
145c4762a1bSJed Brown   PetscScalar    zero = 0.0,one = 1.0;
146c4762a1bSJed Brown   PetscScalar    hx,hy,hxdhy,hydhx;
147c4762a1bSJed Brown   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;
148c4762a1bSJed Brown   PetscScalar    tleft,tright,beta;
149c4762a1bSJed Brown   PetscScalar    **x,**f;
150c4762a1bSJed Brown   Vec            localX;
151c4762a1bSJed Brown   DM             da;
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   PetscFunctionBeginUser;
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localX));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&mx,&my,0,0,0,0,0,0,0,0,0,0));
157c4762a1bSJed Brown   hx    = one/(PetscReal)(mx-1);  hy    = one/(PetscReal)(my-1);
158c4762a1bSJed Brown   hxdhy = hx/hy;               hydhx = hy/hx;
159c4762a1bSJed Brown   tleft = user->tleft;         tright = user->tright;
160c4762a1bSJed Brown   beta  = user->beta;
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /* Get ghost points */
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0));
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,localX,&x));
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,F,&f));
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /* Evaluate function */
170c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
171c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
172c4762a1bSJed Brown       t0 = x[j][i];
173c4762a1bSJed Brown 
174c4762a1bSJed Brown       if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
175c4762a1bSJed Brown 
176c4762a1bSJed Brown         /* general interior volume */
177c4762a1bSJed Brown 
178c4762a1bSJed Brown         tw = x[j][i-1];
179c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
180c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
181c4762a1bSJed Brown         fw = dw*(t0 - tw);
182c4762a1bSJed Brown 
183c4762a1bSJed Brown         te = x[j][i+1];
184c4762a1bSJed Brown         ae = 0.5*(t0 + te);
185c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
186c4762a1bSJed Brown         fe = de*(te - t0);
187c4762a1bSJed Brown 
188c4762a1bSJed Brown         ts = x[j-1][i];
189c4762a1bSJed Brown         as = 0.5*(t0 + ts);
190c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
191c4762a1bSJed Brown         fs = ds*(t0 - ts);
192c4762a1bSJed Brown 
193c4762a1bSJed Brown         tn = x[j+1][i];
194c4762a1bSJed Brown         an = 0.5*(t0 + tn);
195c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
196c4762a1bSJed Brown         fn = dn*(tn - t0);
197c4762a1bSJed Brown 
198c4762a1bSJed Brown       } else if (i == 0) {
199c4762a1bSJed Brown 
200c4762a1bSJed Brown         /* left-hand boundary */
201c4762a1bSJed Brown         tw = tleft;
202c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
203c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
204c4762a1bSJed Brown         fw = dw*(t0 - tw);
205c4762a1bSJed Brown 
206c4762a1bSJed Brown         te = x[j][i+1];
207c4762a1bSJed Brown         ae = 0.5*(t0 + te);
208c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
209c4762a1bSJed Brown         fe = de*(te - t0);
210c4762a1bSJed Brown 
211c4762a1bSJed Brown         if (j > 0) {
212c4762a1bSJed Brown           ts = x[j-1][i];
213c4762a1bSJed Brown           as = 0.5*(t0 + ts);
214c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
215c4762a1bSJed Brown           fs = ds*(t0 - ts);
216c4762a1bSJed Brown         } else fs = zero;
217c4762a1bSJed Brown 
218c4762a1bSJed Brown         if (j < my-1) {
219c4762a1bSJed Brown           tn = x[j+1][i];
220c4762a1bSJed Brown           an = 0.5*(t0 + tn);
221c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
222c4762a1bSJed Brown           fn = dn*(tn - t0);
223c4762a1bSJed Brown         } else fn = zero;
224c4762a1bSJed Brown 
225c4762a1bSJed Brown       } else if (i == mx-1) {
226c4762a1bSJed Brown 
227c4762a1bSJed Brown         /* right-hand boundary */
228c4762a1bSJed Brown         tw = x[j][i-1];
229c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
230c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
231c4762a1bSJed Brown         fw = dw*(t0 - tw);
232c4762a1bSJed Brown 
233c4762a1bSJed Brown         te = tright;
234c4762a1bSJed Brown         ae = 0.5*(t0 + te);
235c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
236c4762a1bSJed Brown         fe = de*(te - t0);
237c4762a1bSJed Brown 
238c4762a1bSJed Brown         if (j > 0) {
239c4762a1bSJed Brown           ts = x[j-1][i];
240c4762a1bSJed Brown           as = 0.5*(t0 + ts);
241c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
242c4762a1bSJed Brown           fs = ds*(t0 - ts);
243c4762a1bSJed Brown         } else fs = zero;
244c4762a1bSJed Brown 
245c4762a1bSJed Brown         if (j < my-1) {
246c4762a1bSJed Brown           tn = x[j+1][i];
247c4762a1bSJed Brown           an = 0.5*(t0 + tn);
248c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
249c4762a1bSJed Brown           fn = dn*(tn - t0);
250c4762a1bSJed Brown         } else fn = zero;
251c4762a1bSJed Brown 
252c4762a1bSJed Brown       } else if (j == 0) {
253c4762a1bSJed Brown 
254c4762a1bSJed Brown         /* bottom boundary,and i <> 0 or mx-1 */
255c4762a1bSJed Brown         tw = x[j][i-1];
256c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
257c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
258c4762a1bSJed Brown         fw = dw*(t0 - tw);
259c4762a1bSJed Brown 
260c4762a1bSJed Brown         te = x[j][i+1];
261c4762a1bSJed Brown         ae = 0.5*(t0 + te);
262c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
263c4762a1bSJed Brown         fe = de*(te - t0);
264c4762a1bSJed Brown 
265c4762a1bSJed Brown         fs = zero;
266c4762a1bSJed Brown 
267c4762a1bSJed Brown         tn = x[j+1][i];
268c4762a1bSJed Brown         an = 0.5*(t0 + tn);
269c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
270c4762a1bSJed Brown         fn = dn*(tn - t0);
271c4762a1bSJed Brown 
272c4762a1bSJed Brown       } else if (j == my-1) {
273c4762a1bSJed Brown 
274c4762a1bSJed Brown         /* top boundary,and i <> 0 or mx-1 */
275c4762a1bSJed Brown         tw = x[j][i-1];
276c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
277c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
278c4762a1bSJed Brown         fw = dw*(t0 - tw);
279c4762a1bSJed Brown 
280c4762a1bSJed Brown         te = x[j][i+1];
281c4762a1bSJed Brown         ae = 0.5*(t0 + te);
282c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
283c4762a1bSJed Brown         fe = de*(te - t0);
284c4762a1bSJed Brown 
285c4762a1bSJed Brown         ts = x[j-1][i];
286c4762a1bSJed Brown         as = 0.5*(t0 + ts);
287c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
288c4762a1bSJed Brown         fs = ds*(t0 - ts);
289c4762a1bSJed Brown 
290c4762a1bSJed Brown         fn = zero;
291c4762a1bSJed Brown 
292c4762a1bSJed Brown       }
293c4762a1bSJed Brown 
294c4762a1bSJed Brown       f[j][i] = -hydhx*(fe-fw) - hxdhy*(fn-fs);
295c4762a1bSJed Brown 
296c4762a1bSJed Brown     }
297c4762a1bSJed Brown   }
2985f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,localX,&x));
2995f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,F,&f));
3005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localX));
3015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops((22.0 + 4.0*POWFLOP)*ym*xm));
302c4762a1bSJed Brown   PetscFunctionReturn(0);
303c4762a1bSJed Brown }
304c4762a1bSJed Brown /* --------------------  Evaluate Jacobian F(x) --------------------- */
305c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X,Mat jac,Mat B,void *ptr)
306c4762a1bSJed Brown {
307c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
308c4762a1bSJed Brown   PetscInt       i,j,mx,my,xs,ys,xm,ym;
309c4762a1bSJed Brown   PetscScalar    one = 1.0,hx,hy,hxdhy,hydhx,t0,tn,ts,te,tw;
310c4762a1bSJed Brown   PetscScalar    dn,ds,de,dw,an,as,ae,aw,bn,bs,be,bw,gn,gs,ge,gw;
311c4762a1bSJed Brown   PetscScalar    tleft,tright,beta,bm1,coef;
312c4762a1bSJed Brown   PetscScalar    v[5],**x;
313c4762a1bSJed Brown   Vec            localX;
314c4762a1bSJed Brown   MatStencil     col[5],row;
315c4762a1bSJed Brown   DM             da;
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   PetscFunctionBeginUser;
3185f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&da));
3195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localX));
3205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,NULL,&mx,&my,0,0,0,0,0,0,0,0,0,0));
321c4762a1bSJed Brown   hx    = one/(PetscReal)(mx-1);  hy     = one/(PetscReal)(my-1);
322c4762a1bSJed Brown   hxdhy = hx/hy;               hydhx  = hy/hx;
323c4762a1bSJed Brown   tleft = user->tleft;         tright = user->tright;
324c4762a1bSJed Brown   beta  = user->beta;          bm1    = user->bm1;          coef = user->coef;
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   /* Get ghost points */
3275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
3285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
3295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0));
3305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,localX,&x));
331c4762a1bSJed Brown 
332c4762a1bSJed Brown   /* Evaluate Jacobian of function */
333c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
334c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
335c4762a1bSJed Brown       t0 = x[j][i];
336c4762a1bSJed Brown 
337c4762a1bSJed Brown       if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
338c4762a1bSJed Brown 
339c4762a1bSJed Brown         /* general interior volume */
340c4762a1bSJed Brown 
341c4762a1bSJed Brown         tw = x[j][i-1];
342c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
343c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
344c4762a1bSJed Brown         /* dw = bw * aw */
345c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
346c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
347c4762a1bSJed Brown 
348c4762a1bSJed Brown         te = x[j][i+1];
349c4762a1bSJed Brown         ae = 0.5*(t0 + te);
350c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
351c4762a1bSJed Brown         /* de = be * ae; */
352c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
353c4762a1bSJed Brown         ge = coef*be*(te - t0);
354c4762a1bSJed Brown 
355c4762a1bSJed Brown         ts = x[j-1][i];
356c4762a1bSJed Brown         as = 0.5*(t0 + ts);
357c4762a1bSJed Brown         bs = PetscPowScalar(as,bm1);
358c4762a1bSJed Brown         /* ds = bs * as; */
359c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
360c4762a1bSJed Brown         gs = coef*bs*(t0 - ts);
361c4762a1bSJed Brown 
362c4762a1bSJed Brown         tn = x[j+1][i];
363c4762a1bSJed Brown         an = 0.5*(t0 + tn);
364c4762a1bSJed Brown         bn = PetscPowScalar(an,bm1);
365c4762a1bSJed Brown         /* dn = bn * an; */
366c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
367c4762a1bSJed Brown         gn = coef*bn*(tn - t0);
368c4762a1bSJed Brown 
369c4762a1bSJed Brown         v[0] = -hxdhy*(ds - gs);                                       col[0].j = j-1;       col[0].i = i;
370c4762a1bSJed Brown         v[1] = -hydhx*(dw - gw);                                       col[1].j = j;         col[1].i = i-1;
371c4762a1bSJed Brown         v[2] = hxdhy*(ds + dn + gs - gn) + hydhx*(dw + de + gw - ge);  col[2].j = row.j = j; col[2].i = row.i = i;
372c4762a1bSJed Brown         v[3] = -hydhx*(de + ge);                                       col[3].j = j;         col[3].i = i+1;
373c4762a1bSJed Brown         v[4] = -hxdhy*(dn + gn);                                       col[4].j = j+1;       col[4].i = i;
3745f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
375c4762a1bSJed Brown 
376c4762a1bSJed Brown       } else if (i == 0) {
377c4762a1bSJed Brown 
378c4762a1bSJed Brown         /* left-hand boundary */
379c4762a1bSJed Brown         tw = tleft;
380c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
381c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
382c4762a1bSJed Brown         /* dw = bw * aw */
383c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
384c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
385c4762a1bSJed Brown 
386c4762a1bSJed Brown         te = x[j][i + 1];
387c4762a1bSJed Brown         ae = 0.5*(t0 + te);
388c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
389c4762a1bSJed Brown         /* de = be * ae; */
390c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
391c4762a1bSJed Brown         ge = coef*be*(te - t0);
392c4762a1bSJed Brown 
393c4762a1bSJed Brown         /* left-hand bottom boundary */
394c4762a1bSJed Brown         if (j == 0) {
395c4762a1bSJed Brown 
396c4762a1bSJed Brown           tn = x[j+1][i];
397c4762a1bSJed Brown           an = 0.5*(t0 + tn);
398c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
399c4762a1bSJed Brown           /* dn = bn * an; */
400c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
401c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
402c4762a1bSJed Brown 
403c4762a1bSJed Brown           v[0] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[0].j = row.j = j; col[0].i = row.i = i;
404c4762a1bSJed Brown           v[1] = -hydhx*(de + ge);                            col[1].j = j;         col[1].i = i+1;
405c4762a1bSJed Brown           v[2] = -hxdhy*(dn + gn);                            col[2].j = j+1;       col[2].i = i;
4065f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
407c4762a1bSJed Brown 
408c4762a1bSJed Brown           /* left-hand interior boundary */
409c4762a1bSJed Brown         } else if (j < my-1) {
410c4762a1bSJed Brown 
411c4762a1bSJed Brown           ts = x[j-1][i];
412c4762a1bSJed Brown           as = 0.5*(t0 + ts);
413c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
414c4762a1bSJed Brown           /* ds = bs * as; */
415c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
416c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
417c4762a1bSJed Brown 
418c4762a1bSJed Brown           tn = x[j+1][i];
419c4762a1bSJed Brown           an = 0.5*(t0 + tn);
420c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
421c4762a1bSJed Brown           /* dn = bn * an; */
422c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
423c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
424c4762a1bSJed Brown 
425c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                                       col[0].j = j-1;       col[0].i = i;
426c4762a1bSJed Brown           v[1] = hxdhy*(ds + dn + gs - gn) + hydhx*(dw + de + gw - ge);  col[1].j = row.j = j; col[1].i = row.i = i;
427c4762a1bSJed Brown           v[2] = -hydhx*(de + ge);                                       col[2].j = j;         col[2].i = i+1;
428c4762a1bSJed Brown           v[3] = -hxdhy*(dn + gn);                                       col[3].j = j+1;       col[3].i = i;
4295f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
430c4762a1bSJed Brown           /* left-hand top boundary */
431c4762a1bSJed Brown         } else {
432c4762a1bSJed Brown 
433c4762a1bSJed Brown           ts = x[j-1][i];
434c4762a1bSJed Brown           as = 0.5*(t0 + ts);
435c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
436c4762a1bSJed Brown           /* ds = bs * as; */
437c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
438c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
439c4762a1bSJed Brown 
440c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                             col[0].j = j-1;       col[0].i = i;
441c4762a1bSJed Brown           v[1] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge);  col[1].j = row.j = j; col[1].i = row.i = i;
442c4762a1bSJed Brown           v[2] = -hydhx*(de + ge);                             col[2].j = j;         col[2].i = i+1;
4435f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
444c4762a1bSJed Brown         }
445c4762a1bSJed Brown 
446c4762a1bSJed Brown       } else if (i == mx-1) {
447c4762a1bSJed Brown 
448c4762a1bSJed Brown         /* right-hand boundary */
449c4762a1bSJed Brown         tw = x[j][i-1];
450c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
451c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
452c4762a1bSJed Brown         /* dw = bw * aw */
453c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
454c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
455c4762a1bSJed Brown 
456c4762a1bSJed Brown         te = tright;
457c4762a1bSJed Brown         ae = 0.5*(t0 + te);
458c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
459c4762a1bSJed Brown         /* de = be * ae; */
460c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
461c4762a1bSJed Brown         ge = coef*be*(te - t0);
462c4762a1bSJed Brown 
463c4762a1bSJed Brown         /* right-hand bottom boundary */
464c4762a1bSJed Brown         if (j == 0) {
465c4762a1bSJed Brown 
466c4762a1bSJed Brown           tn = x[j+1][i];
467c4762a1bSJed Brown           an = 0.5*(t0 + tn);
468c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
469c4762a1bSJed Brown           /* dn = bn * an; */
470c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
471c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
472c4762a1bSJed Brown 
473c4762a1bSJed Brown           v[0] = -hydhx*(dw - gw);                            col[0].j = j;         col[0].i = i-1;
474c4762a1bSJed Brown           v[1] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[1].j = row.j = j; col[1].i = row.i = i;
475c4762a1bSJed Brown           v[2] = -hxdhy*(dn + gn);                            col[2].j = j+1;       col[2].i = i;
4765f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
477c4762a1bSJed Brown 
478c4762a1bSJed Brown           /* right-hand interior boundary */
479c4762a1bSJed Brown         } else if (j < my-1) {
480c4762a1bSJed Brown 
481c4762a1bSJed Brown           ts = x[j-1][i];
482c4762a1bSJed Brown           as = 0.5*(t0 + ts);
483c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
484c4762a1bSJed Brown           /* ds = bs * as; */
485c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
486c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
487c4762a1bSJed Brown 
488c4762a1bSJed Brown           tn = x[j+1][i];
489c4762a1bSJed Brown           an = 0.5*(t0 + tn);
490c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
491c4762a1bSJed Brown           /* dn = bn * an; */
492c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
493c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
494c4762a1bSJed Brown 
495c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                                       col[0].j = j-1;       col[0].i = i;
496c4762a1bSJed Brown           v[1] = -hydhx*(dw - gw);                                       col[1].j = j;         col[1].i = i-1;
497c4762a1bSJed Brown           v[2] = hxdhy*(ds + dn + gs - gn) + hydhx*(dw + de + gw - ge);  col[2].j = row.j = j; col[2].i = row.i = i;
498c4762a1bSJed Brown           v[3] = -hxdhy*(dn + gn);                                       col[3].j = j+1;       col[3].i = i;
4995f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
500c4762a1bSJed Brown         /* right-hand top boundary */
501c4762a1bSJed Brown         } else {
502c4762a1bSJed Brown 
503c4762a1bSJed Brown           ts = x[j-1][i];
504c4762a1bSJed Brown           as = 0.5*(t0 + ts);
505c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
506c4762a1bSJed Brown           /* ds = bs * as; */
507c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
508c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
509c4762a1bSJed Brown 
510c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                             col[0].j = j-1;       col[0].i = i;
511c4762a1bSJed Brown           v[1] = -hydhx*(dw - gw);                             col[1].j = j;         col[1].i = i-1;
512c4762a1bSJed Brown           v[2] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge);  col[2].j = row.j = j; col[2].i = row.i = i;
5135f80ce2aSJacob Faibussowitsch           CHKERRQ(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
514c4762a1bSJed Brown         }
515c4762a1bSJed Brown 
516c4762a1bSJed Brown         /* bottom boundary,and i <> 0 or mx-1 */
517c4762a1bSJed Brown       } else if (j == 0) {
518c4762a1bSJed Brown 
519c4762a1bSJed Brown         tw = x[j][i-1];
520c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
521c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
522c4762a1bSJed Brown         /* dw = bw * aw */
523c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
524c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
525c4762a1bSJed Brown 
526c4762a1bSJed Brown         te = x[j][i+1];
527c4762a1bSJed Brown         ae = 0.5*(t0 + te);
528c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
529c4762a1bSJed Brown         /* de = be * ae; */
530c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
531c4762a1bSJed Brown         ge = coef*be*(te - t0);
532c4762a1bSJed Brown 
533c4762a1bSJed Brown         tn = x[j+1][i];
534c4762a1bSJed Brown         an = 0.5*(t0 + tn);
535c4762a1bSJed Brown         bn = PetscPowScalar(an,bm1);
536c4762a1bSJed Brown         /* dn = bn * an; */
537c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
538c4762a1bSJed Brown         gn = coef*bn*(tn - t0);
539c4762a1bSJed Brown 
540c4762a1bSJed Brown         v[0] = -hydhx*(dw - gw);                            col[0].j = j;         col[0].i = i-1;
541c4762a1bSJed Brown         v[1] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[1].j = row.j = j; col[1].i = row.i = i;
542c4762a1bSJed Brown         v[2] = -hydhx*(de + ge);                            col[2].j = j;         col[2].i = i+1;
543c4762a1bSJed Brown         v[3] = -hxdhy*(dn + gn);                            col[3].j = j+1;       col[3].i = i;
5445f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
545c4762a1bSJed Brown 
546c4762a1bSJed Brown         /* top boundary,and i <> 0 or mx-1 */
547c4762a1bSJed Brown       } else if (j == my-1) {
548c4762a1bSJed Brown 
549c4762a1bSJed Brown         tw = x[j][i-1];
550c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
551c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
552c4762a1bSJed Brown         /* dw = bw * aw */
553c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
554c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
555c4762a1bSJed Brown 
556c4762a1bSJed Brown         te = x[j][i+1];
557c4762a1bSJed Brown         ae = 0.5*(t0 + te);
558c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
559c4762a1bSJed Brown         /* de = be * ae; */
560c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
561c4762a1bSJed Brown         ge = coef*be*(te - t0);
562c4762a1bSJed Brown 
563c4762a1bSJed Brown         ts = x[j-1][i];
564c4762a1bSJed Brown         as = 0.5*(t0 + ts);
565c4762a1bSJed Brown         bs = PetscPowScalar(as,bm1);
566c4762a1bSJed Brown         /* ds = bs * as; */
567c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
568c4762a1bSJed Brown         gs = coef*bs*(t0 - ts);
569c4762a1bSJed Brown 
570c4762a1bSJed Brown         v[0] = -hxdhy*(ds - gs);                             col[0].j = j-1;       col[0].i = i;
571c4762a1bSJed Brown         v[1] = -hydhx*(dw - gw);                             col[1].j = j;         col[1].i = i-1;
572c4762a1bSJed Brown         v[2] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge);  col[2].j = row.j = j; col[2].i = row.i = i;
573c4762a1bSJed Brown         v[3] = -hydhx*(de + ge);                             col[3].j = j;         col[3].i = i+1;
5745f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
575c4762a1bSJed Brown 
576c4762a1bSJed Brown       }
577c4762a1bSJed Brown     }
578c4762a1bSJed Brown   }
5795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
5805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,localX,&x));
5815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
5825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localX));
583c4762a1bSJed Brown   if (jac != B) {
5845f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
5855f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
586c4762a1bSJed Brown   }
587c4762a1bSJed Brown 
5885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops((41.0 + 8.0*POWFLOP)*xm*ym));
589c4762a1bSJed Brown   PetscFunctionReturn(0);
590c4762a1bSJed Brown }
591c4762a1bSJed Brown 
592c4762a1bSJed Brown /*TEST
593c4762a1bSJed Brown 
594c4762a1bSJed Brown    test:
59541ba4c6cSHeeho Park       suffix: 1
596c4762a1bSJed Brown       args: -pc_type mg -ksp_type fgmres -da_refine 2 -pc_mg_galerkin pmat -snes_view
597c4762a1bSJed Brown       requires: !single
598c4762a1bSJed Brown 
59941ba4c6cSHeeho Park    test:
60041ba4c6cSHeeho Park       suffix: 2
60141ba4c6cSHeeho Park       args: -pc_type mg -ksp_type fgmres -da_refine 2 -pc_mg_galerkin pmat -snes_view -snes_type newtontrdc
60241ba4c6cSHeeho Park       requires: !single
60341ba4c6cSHeeho Park 
604c4762a1bSJed Brown TEST*/
605