xref: /petsc/src/snes/tutorials/ex18.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 /*
15c4762a1bSJed Brown 
16c4762a1bSJed Brown     This example models the partial differential equation
17c4762a1bSJed Brown 
18c4762a1bSJed Brown          - Div(alpha* T^beta (GRAD T)) = 0.
19c4762a1bSJed Brown 
20c4762a1bSJed Brown     where beta = 2.5 and alpha = 1.0
21c4762a1bSJed Brown 
22c4762a1bSJed Brown     BC: T_left = 1.0, T_right = 0.1, dT/dn_top = dTdn_bottom = 0.
23c4762a1bSJed Brown 
24c4762a1bSJed Brown     in the unit square, which is uniformly discretized in each of x and
25c4762a1bSJed Brown     y in this simple encoding.  The degrees of freedom are cell centered.
26c4762a1bSJed Brown 
27c4762a1bSJed Brown     A finite volume approximation with the usual 5-point stencil
28c4762a1bSJed Brown     is used to discretize the boundary value problem to obtain a
29c4762a1bSJed Brown     nonlinear system of equations.
30c4762a1bSJed Brown 
31c4762a1bSJed Brown     This code was contributed by David Keyes
32c4762a1bSJed Brown 
33c4762a1bSJed Brown */
34c4762a1bSJed Brown 
35c4762a1bSJed Brown #include <petscsnes.h>
36c4762a1bSJed Brown #include <petscdm.h>
37c4762a1bSJed Brown #include <petscdmda.h>
38c4762a1bSJed Brown 
39c4762a1bSJed Brown /* User-defined application context */
40c4762a1bSJed Brown 
41c4762a1bSJed Brown typedef struct {
42c4762a1bSJed Brown   PetscReal tleft,tright;    /* Dirichlet boundary conditions */
43c4762a1bSJed Brown   PetscReal beta,bm1,coef;   /* nonlinear diffusivity parameterizations */
44c4762a1bSJed Brown } AppCtx;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown #define POWFLOP 5 /* assume a pow() takes five flops */
47c4762a1bSJed Brown 
48c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(SNES,Vec,void*);
49c4762a1bSJed Brown extern PetscErrorCode FormFunction(SNES,Vec,Vec,void*);
50c4762a1bSJed Brown extern PetscErrorCode FormJacobian(SNES,Vec,Mat,Mat,void*);
51c4762a1bSJed Brown 
52c4762a1bSJed Brown int main(int argc,char **argv)
53c4762a1bSJed Brown {
54c4762a1bSJed Brown   SNES           snes;
55c4762a1bSJed Brown   AppCtx         user;
56c4762a1bSJed Brown   PetscInt       its,lits;
57c4762a1bSJed Brown   PetscReal      litspit;
58c4762a1bSJed Brown   DM             da;
59c4762a1bSJed Brown 
60*327415f7SBarry Smith   PetscFunctionBeginUser;
619566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   /* set problem parameters */
64c4762a1bSJed Brown   user.tleft  = 1.0;
65c4762a1bSJed Brown   user.tright = 0.1;
66c4762a1bSJed Brown   user.beta   = 2.5;
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-tleft",&user.tleft,NULL));
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-tright",&user.tright,NULL));
699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-beta",&user.beta,NULL));
70c4762a1bSJed Brown   user.bm1    = user.beta - 1.0;
71c4762a1bSJed Brown   user.coef   = user.beta/2.0;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /*
74c4762a1bSJed Brown       Create the multilevel DM data structure
75c4762a1bSJed Brown   */
769566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /*
79c4762a1bSJed Brown       Set the DMDA (grid structure) for the grids.
80c4762a1bSJed Brown   */
819566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,5,5,PETSC_DECIDE,PETSC_DECIDE,1,1,0,0,&da));
829566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
839566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
849566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da,&user));
859566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes,(DM)da));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /*
88c4762a1bSJed Brown      Create the nonlinear solver, and tell it the functions to use
89c4762a1bSJed Brown   */
909566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes,NULL,FormFunction,&user));
919566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes,NULL,NULL,FormJacobian,&user));
929566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
939566063dSJacob Faibussowitsch   PetscCall(SNESSetComputeInitialGuess(snes,FormInitialGuess,NULL));
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,NULL,NULL));
969566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&its));
979566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(snes,&lits));
98c4762a1bSJed Brown   litspit = ((PetscReal)lits)/((PetscReal)its);
9963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %" PetscInt_FMT "\n",its));
10063a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of Linear iterations = %" PetscInt_FMT "\n",lits));
1019566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Average Linear its / SNES = %e\n",(double)litspit));
102c4762a1bSJed Brown 
1039566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1049566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1059566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
106b122ec5aSJacob Faibussowitsch   return 0;
107c4762a1bSJed Brown }
108c4762a1bSJed Brown /* --------------------  Form initial approximation ----------------- */
109c4762a1bSJed Brown PetscErrorCode FormInitialGuess(SNES snes,Vec X,void *ctx)
110c4762a1bSJed Brown {
111c4762a1bSJed Brown   AppCtx         *user;
112c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym;
113c4762a1bSJed Brown   PetscReal      tleft;
114c4762a1bSJed Brown   PetscScalar    **x;
115c4762a1bSJed Brown   DM             da;
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   PetscFunctionBeginUser;
1189566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&da));
1199566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da,&user));
120c4762a1bSJed Brown   tleft = user->tleft;
121c4762a1bSJed Brown   /* Get ghost points */
1229566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0));
1239566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,X,&x));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Compute initial guess */
126c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
127c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
128c4762a1bSJed Brown       x[j][i] = tleft;
129c4762a1bSJed Brown     }
130c4762a1bSJed Brown   }
1319566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,X,&x));
132c4762a1bSJed Brown   PetscFunctionReturn(0);
133c4762a1bSJed Brown }
134c4762a1bSJed Brown /* --------------------  Evaluate Function F(x) --------------------- */
135c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
136c4762a1bSJed Brown {
137c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
138c4762a1bSJed Brown   PetscInt       i,j,mx,my,xs,ys,xm,ym;
139c4762a1bSJed Brown   PetscScalar    zero = 0.0,one = 1.0;
140c4762a1bSJed Brown   PetscScalar    hx,hy,hxdhy,hydhx;
141c4762a1bSJed 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;
142c4762a1bSJed Brown   PetscScalar    tleft,tright,beta;
143c4762a1bSJed Brown   PetscScalar    **x,**f;
144c4762a1bSJed Brown   Vec            localX;
145c4762a1bSJed Brown   DM             da;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   PetscFunctionBeginUser;
1489566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&da));
1499566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localX));
1509566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,NULL,&mx,&my,0,0,0,0,0,0,0,0,0,0));
151c4762a1bSJed Brown   hx    = one/(PetscReal)(mx-1);  hy    = one/(PetscReal)(my-1);
152c4762a1bSJed Brown   hxdhy = hx/hy;               hydhx = hy/hx;
153c4762a1bSJed Brown   tleft = user->tleft;         tright = user->tright;
154c4762a1bSJed Brown   beta  = user->beta;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /* Get ghost points */
1579566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
1589566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
1599566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0));
1609566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,localX,&x));
1619566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,F,&f));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* Evaluate function */
164c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
165c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
166c4762a1bSJed Brown       t0 = x[j][i];
167c4762a1bSJed Brown 
168c4762a1bSJed Brown       if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
169c4762a1bSJed Brown 
170c4762a1bSJed Brown         /* general interior volume */
171c4762a1bSJed Brown 
172c4762a1bSJed Brown         tw = x[j][i-1];
173c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
174c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
175c4762a1bSJed Brown         fw = dw*(t0 - tw);
176c4762a1bSJed Brown 
177c4762a1bSJed Brown         te = x[j][i+1];
178c4762a1bSJed Brown         ae = 0.5*(t0 + te);
179c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
180c4762a1bSJed Brown         fe = de*(te - t0);
181c4762a1bSJed Brown 
182c4762a1bSJed Brown         ts = x[j-1][i];
183c4762a1bSJed Brown         as = 0.5*(t0 + ts);
184c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
185c4762a1bSJed Brown         fs = ds*(t0 - ts);
186c4762a1bSJed Brown 
187c4762a1bSJed Brown         tn = x[j+1][i];
188c4762a1bSJed Brown         an = 0.5*(t0 + tn);
189c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
190c4762a1bSJed Brown         fn = dn*(tn - t0);
191c4762a1bSJed Brown 
192c4762a1bSJed Brown       } else if (i == 0) {
193c4762a1bSJed Brown 
194c4762a1bSJed Brown         /* left-hand boundary */
195c4762a1bSJed Brown         tw = tleft;
196c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
197c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
198c4762a1bSJed Brown         fw = dw*(t0 - tw);
199c4762a1bSJed Brown 
200c4762a1bSJed Brown         te = x[j][i+1];
201c4762a1bSJed Brown         ae = 0.5*(t0 + te);
202c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
203c4762a1bSJed Brown         fe = de*(te - t0);
204c4762a1bSJed Brown 
205c4762a1bSJed Brown         if (j > 0) {
206c4762a1bSJed Brown           ts = x[j-1][i];
207c4762a1bSJed Brown           as = 0.5*(t0 + ts);
208c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
209c4762a1bSJed Brown           fs = ds*(t0 - ts);
210c4762a1bSJed Brown         } else fs = zero;
211c4762a1bSJed Brown 
212c4762a1bSJed Brown         if (j < my-1) {
213c4762a1bSJed Brown           tn = x[j+1][i];
214c4762a1bSJed Brown           an = 0.5*(t0 + tn);
215c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
216c4762a1bSJed Brown           fn = dn*(tn - t0);
217c4762a1bSJed Brown         } else fn = zero;
218c4762a1bSJed Brown 
219c4762a1bSJed Brown       } else if (i == mx-1) {
220c4762a1bSJed Brown 
221c4762a1bSJed Brown         /* right-hand boundary */
222c4762a1bSJed Brown         tw = x[j][i-1];
223c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
224c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
225c4762a1bSJed Brown         fw = dw*(t0 - tw);
226c4762a1bSJed Brown 
227c4762a1bSJed Brown         te = tright;
228c4762a1bSJed Brown         ae = 0.5*(t0 + te);
229c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
230c4762a1bSJed Brown         fe = de*(te - t0);
231c4762a1bSJed Brown 
232c4762a1bSJed Brown         if (j > 0) {
233c4762a1bSJed Brown           ts = x[j-1][i];
234c4762a1bSJed Brown           as = 0.5*(t0 + ts);
235c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
236c4762a1bSJed Brown           fs = ds*(t0 - ts);
237c4762a1bSJed Brown         } else fs = zero;
238c4762a1bSJed Brown 
239c4762a1bSJed Brown         if (j < my-1) {
240c4762a1bSJed Brown           tn = x[j+1][i];
241c4762a1bSJed Brown           an = 0.5*(t0 + tn);
242c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
243c4762a1bSJed Brown           fn = dn*(tn - t0);
244c4762a1bSJed Brown         } else fn = zero;
245c4762a1bSJed Brown 
246c4762a1bSJed Brown       } else if (j == 0) {
247c4762a1bSJed Brown 
248c4762a1bSJed Brown         /* bottom boundary,and i <> 0 or mx-1 */
249c4762a1bSJed Brown         tw = x[j][i-1];
250c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
251c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
252c4762a1bSJed Brown         fw = dw*(t0 - tw);
253c4762a1bSJed Brown 
254c4762a1bSJed Brown         te = x[j][i+1];
255c4762a1bSJed Brown         ae = 0.5*(t0 + te);
256c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
257c4762a1bSJed Brown         fe = de*(te - t0);
258c4762a1bSJed Brown 
259c4762a1bSJed Brown         fs = zero;
260c4762a1bSJed Brown 
261c4762a1bSJed Brown         tn = x[j+1][i];
262c4762a1bSJed Brown         an = 0.5*(t0 + tn);
263c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
264c4762a1bSJed Brown         fn = dn*(tn - t0);
265c4762a1bSJed Brown 
266c4762a1bSJed Brown       } else if (j == my-1) {
267c4762a1bSJed Brown 
268c4762a1bSJed Brown         /* top boundary,and i <> 0 or mx-1 */
269c4762a1bSJed Brown         tw = x[j][i-1];
270c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
271c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
272c4762a1bSJed Brown         fw = dw*(t0 - tw);
273c4762a1bSJed Brown 
274c4762a1bSJed Brown         te = x[j][i+1];
275c4762a1bSJed Brown         ae = 0.5*(t0 + te);
276c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
277c4762a1bSJed Brown         fe = de*(te - t0);
278c4762a1bSJed Brown 
279c4762a1bSJed Brown         ts = x[j-1][i];
280c4762a1bSJed Brown         as = 0.5*(t0 + ts);
281c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
282c4762a1bSJed Brown         fs = ds*(t0 - ts);
283c4762a1bSJed Brown 
284c4762a1bSJed Brown         fn = zero;
285c4762a1bSJed Brown 
286c4762a1bSJed Brown       }
287c4762a1bSJed Brown 
288c4762a1bSJed Brown       f[j][i] = -hydhx*(fe-fw) - hxdhy*(fn-fs);
289c4762a1bSJed Brown 
290c4762a1bSJed Brown     }
291c4762a1bSJed Brown   }
2929566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,localX,&x));
2939566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,F,&f));
2949566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localX));
2959566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops((22.0 + 4.0*POWFLOP)*ym*xm));
296c4762a1bSJed Brown   PetscFunctionReturn(0);
297c4762a1bSJed Brown }
298c4762a1bSJed Brown /* --------------------  Evaluate Jacobian F(x) --------------------- */
299c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X,Mat jac,Mat B,void *ptr)
300c4762a1bSJed Brown {
301c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
302c4762a1bSJed Brown   PetscInt       i,j,mx,my,xs,ys,xm,ym;
303c4762a1bSJed Brown   PetscScalar    one = 1.0,hx,hy,hxdhy,hydhx,t0,tn,ts,te,tw;
304c4762a1bSJed Brown   PetscScalar    dn,ds,de,dw,an,as,ae,aw,bn,bs,be,bw,gn,gs,ge,gw;
305c4762a1bSJed Brown   PetscScalar    tleft,tright,beta,bm1,coef;
306c4762a1bSJed Brown   PetscScalar    v[5],**x;
307c4762a1bSJed Brown   Vec            localX;
308c4762a1bSJed Brown   MatStencil     col[5],row;
309c4762a1bSJed Brown   DM             da;
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   PetscFunctionBeginUser;
3129566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&da));
3139566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localX));
3149566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,NULL,&mx,&my,0,0,0,0,0,0,0,0,0,0));
315c4762a1bSJed Brown   hx    = one/(PetscReal)(mx-1);  hy     = one/(PetscReal)(my-1);
316c4762a1bSJed Brown   hxdhy = hx/hy;               hydhx  = hy/hx;
317c4762a1bSJed Brown   tleft = user->tleft;         tright = user->tright;
318c4762a1bSJed Brown   beta  = user->beta;          bm1    = user->bm1;          coef = user->coef;
319c4762a1bSJed Brown 
320c4762a1bSJed Brown   /* Get ghost points */
3219566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
3229566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
3239566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0));
3249566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,localX,&x));
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   /* Evaluate Jacobian of function */
327c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
328c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
329c4762a1bSJed Brown       t0 = x[j][i];
330c4762a1bSJed Brown 
331c4762a1bSJed Brown       if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
332c4762a1bSJed Brown 
333c4762a1bSJed Brown         /* general interior volume */
334c4762a1bSJed Brown 
335c4762a1bSJed Brown         tw = x[j][i-1];
336c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
337c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
338c4762a1bSJed Brown         /* dw = bw * aw */
339c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
340c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
341c4762a1bSJed Brown 
342c4762a1bSJed Brown         te = x[j][i+1];
343c4762a1bSJed Brown         ae = 0.5*(t0 + te);
344c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
345c4762a1bSJed Brown         /* de = be * ae; */
346c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
347c4762a1bSJed Brown         ge = coef*be*(te - t0);
348c4762a1bSJed Brown 
349c4762a1bSJed Brown         ts = x[j-1][i];
350c4762a1bSJed Brown         as = 0.5*(t0 + ts);
351c4762a1bSJed Brown         bs = PetscPowScalar(as,bm1);
352c4762a1bSJed Brown         /* ds = bs * as; */
353c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
354c4762a1bSJed Brown         gs = coef*bs*(t0 - ts);
355c4762a1bSJed Brown 
356c4762a1bSJed Brown         tn = x[j+1][i];
357c4762a1bSJed Brown         an = 0.5*(t0 + tn);
358c4762a1bSJed Brown         bn = PetscPowScalar(an,bm1);
359c4762a1bSJed Brown         /* dn = bn * an; */
360c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
361c4762a1bSJed Brown         gn = coef*bn*(tn - t0);
362c4762a1bSJed Brown 
363c4762a1bSJed Brown         v[0] = -hxdhy*(ds - gs);                                       col[0].j = j-1;       col[0].i = i;
364c4762a1bSJed Brown         v[1] = -hydhx*(dw - gw);                                       col[1].j = j;         col[1].i = i-1;
365c4762a1bSJed 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;
366c4762a1bSJed Brown         v[3] = -hydhx*(de + ge);                                       col[3].j = j;         col[3].i = i+1;
367c4762a1bSJed Brown         v[4] = -hxdhy*(dn + gn);                                       col[4].j = j+1;       col[4].i = i;
3689566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
369c4762a1bSJed Brown 
370c4762a1bSJed Brown       } else if (i == 0) {
371c4762a1bSJed Brown 
372c4762a1bSJed Brown         /* left-hand boundary */
373c4762a1bSJed Brown         tw = tleft;
374c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
375c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
376c4762a1bSJed Brown         /* dw = bw * aw */
377c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
378c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
379c4762a1bSJed Brown 
380c4762a1bSJed Brown         te = x[j][i + 1];
381c4762a1bSJed Brown         ae = 0.5*(t0 + te);
382c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
383c4762a1bSJed Brown         /* de = be * ae; */
384c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
385c4762a1bSJed Brown         ge = coef*be*(te - t0);
386c4762a1bSJed Brown 
387c4762a1bSJed Brown         /* left-hand bottom boundary */
388c4762a1bSJed Brown         if (j == 0) {
389c4762a1bSJed Brown 
390c4762a1bSJed Brown           tn = x[j+1][i];
391c4762a1bSJed Brown           an = 0.5*(t0 + tn);
392c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
393c4762a1bSJed Brown           /* dn = bn * an; */
394c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
395c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
396c4762a1bSJed Brown 
397c4762a1bSJed Brown           v[0] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[0].j = row.j = j; col[0].i = row.i = i;
398c4762a1bSJed Brown           v[1] = -hydhx*(de + ge);                            col[1].j = j;         col[1].i = i+1;
399c4762a1bSJed Brown           v[2] = -hxdhy*(dn + gn);                            col[2].j = j+1;       col[2].i = i;
4009566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
401c4762a1bSJed Brown 
402c4762a1bSJed Brown           /* left-hand interior boundary */
403c4762a1bSJed Brown         } else if (j < my-1) {
404c4762a1bSJed Brown 
405c4762a1bSJed Brown           ts = x[j-1][i];
406c4762a1bSJed Brown           as = 0.5*(t0 + ts);
407c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
408c4762a1bSJed Brown           /* ds = bs * as; */
409c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
410c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
411c4762a1bSJed Brown 
412c4762a1bSJed Brown           tn = x[j+1][i];
413c4762a1bSJed Brown           an = 0.5*(t0 + tn);
414c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
415c4762a1bSJed Brown           /* dn = bn * an; */
416c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
417c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
418c4762a1bSJed Brown 
419c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                                       col[0].j = j-1;       col[0].i = i;
420c4762a1bSJed 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;
421c4762a1bSJed Brown           v[2] = -hydhx*(de + ge);                                       col[2].j = j;         col[2].i = i+1;
422c4762a1bSJed Brown           v[3] = -hxdhy*(dn + gn);                                       col[3].j = j+1;       col[3].i = i;
4239566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
424c4762a1bSJed Brown           /* left-hand top boundary */
425c4762a1bSJed Brown         } else {
426c4762a1bSJed Brown 
427c4762a1bSJed Brown           ts = x[j-1][i];
428c4762a1bSJed Brown           as = 0.5*(t0 + ts);
429c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
430c4762a1bSJed Brown           /* ds = bs * as; */
431c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
432c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
433c4762a1bSJed Brown 
434c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                             col[0].j = j-1;       col[0].i = i;
435c4762a1bSJed Brown           v[1] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge);  col[1].j = row.j = j; col[1].i = row.i = i;
436c4762a1bSJed Brown           v[2] = -hydhx*(de + ge);                             col[2].j = j;         col[2].i = i+1;
4379566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
438c4762a1bSJed Brown         }
439c4762a1bSJed Brown 
440c4762a1bSJed Brown       } else if (i == mx-1) {
441c4762a1bSJed Brown 
442c4762a1bSJed Brown         /* right-hand boundary */
443c4762a1bSJed Brown         tw = x[j][i-1];
444c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
445c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
446c4762a1bSJed Brown         /* dw = bw * aw */
447c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
448c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
449c4762a1bSJed Brown 
450c4762a1bSJed Brown         te = tright;
451c4762a1bSJed Brown         ae = 0.5*(t0 + te);
452c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
453c4762a1bSJed Brown         /* de = be * ae; */
454c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
455c4762a1bSJed Brown         ge = coef*be*(te - t0);
456c4762a1bSJed Brown 
457c4762a1bSJed Brown         /* right-hand bottom boundary */
458c4762a1bSJed Brown         if (j == 0) {
459c4762a1bSJed Brown 
460c4762a1bSJed Brown           tn = x[j+1][i];
461c4762a1bSJed Brown           an = 0.5*(t0 + tn);
462c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
463c4762a1bSJed Brown           /* dn = bn * an; */
464c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
465c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
466c4762a1bSJed Brown 
467c4762a1bSJed Brown           v[0] = -hydhx*(dw - gw);                            col[0].j = j;         col[0].i = i-1;
468c4762a1bSJed Brown           v[1] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[1].j = row.j = j; col[1].i = row.i = i;
469c4762a1bSJed Brown           v[2] = -hxdhy*(dn + gn);                            col[2].j = j+1;       col[2].i = i;
4709566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
471c4762a1bSJed Brown 
472c4762a1bSJed Brown           /* right-hand interior boundary */
473c4762a1bSJed Brown         } else if (j < my-1) {
474c4762a1bSJed Brown 
475c4762a1bSJed Brown           ts = x[j-1][i];
476c4762a1bSJed Brown           as = 0.5*(t0 + ts);
477c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
478c4762a1bSJed Brown           /* ds = bs * as; */
479c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
480c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
481c4762a1bSJed Brown 
482c4762a1bSJed Brown           tn = x[j+1][i];
483c4762a1bSJed Brown           an = 0.5*(t0 + tn);
484c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
485c4762a1bSJed Brown           /* dn = bn * an; */
486c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
487c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
488c4762a1bSJed Brown 
489c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                                       col[0].j = j-1;       col[0].i = i;
490c4762a1bSJed Brown           v[1] = -hydhx*(dw - gw);                                       col[1].j = j;         col[1].i = i-1;
491c4762a1bSJed 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;
492c4762a1bSJed Brown           v[3] = -hxdhy*(dn + gn);                                       col[3].j = j+1;       col[3].i = i;
4939566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
494c4762a1bSJed Brown         /* right-hand top boundary */
495c4762a1bSJed Brown         } else {
496c4762a1bSJed Brown 
497c4762a1bSJed Brown           ts = x[j-1][i];
498c4762a1bSJed Brown           as = 0.5*(t0 + ts);
499c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
500c4762a1bSJed Brown           /* ds = bs * as; */
501c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
502c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
503c4762a1bSJed Brown 
504c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                             col[0].j = j-1;       col[0].i = i;
505c4762a1bSJed Brown           v[1] = -hydhx*(dw - gw);                             col[1].j = j;         col[1].i = i-1;
506c4762a1bSJed Brown           v[2] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge);  col[2].j = row.j = j; col[2].i = row.i = i;
5079566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
508c4762a1bSJed Brown         }
509c4762a1bSJed Brown 
510c4762a1bSJed Brown         /* bottom boundary,and i <> 0 or mx-1 */
511c4762a1bSJed Brown       } else if (j == 0) {
512c4762a1bSJed Brown 
513c4762a1bSJed Brown         tw = x[j][i-1];
514c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
515c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
516c4762a1bSJed Brown         /* dw = bw * aw */
517c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
518c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
519c4762a1bSJed Brown 
520c4762a1bSJed Brown         te = x[j][i+1];
521c4762a1bSJed Brown         ae = 0.5*(t0 + te);
522c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
523c4762a1bSJed Brown         /* de = be * ae; */
524c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
525c4762a1bSJed Brown         ge = coef*be*(te - t0);
526c4762a1bSJed Brown 
527c4762a1bSJed Brown         tn = x[j+1][i];
528c4762a1bSJed Brown         an = 0.5*(t0 + tn);
529c4762a1bSJed Brown         bn = PetscPowScalar(an,bm1);
530c4762a1bSJed Brown         /* dn = bn * an; */
531c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
532c4762a1bSJed Brown         gn = coef*bn*(tn - t0);
533c4762a1bSJed Brown 
534c4762a1bSJed Brown         v[0] = -hydhx*(dw - gw);                            col[0].j = j;         col[0].i = i-1;
535c4762a1bSJed Brown         v[1] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[1].j = row.j = j; col[1].i = row.i = i;
536c4762a1bSJed Brown         v[2] = -hydhx*(de + ge);                            col[2].j = j;         col[2].i = i+1;
537c4762a1bSJed Brown         v[3] = -hxdhy*(dn + gn);                            col[3].j = j+1;       col[3].i = i;
5389566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
539c4762a1bSJed Brown 
540c4762a1bSJed Brown         /* top boundary,and i <> 0 or mx-1 */
541c4762a1bSJed Brown       } else if (j == my-1) {
542c4762a1bSJed Brown 
543c4762a1bSJed Brown         tw = x[j][i-1];
544c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
545c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
546c4762a1bSJed Brown         /* dw = bw * aw */
547c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
548c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
549c4762a1bSJed Brown 
550c4762a1bSJed Brown         te = x[j][i+1];
551c4762a1bSJed Brown         ae = 0.5*(t0 + te);
552c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
553c4762a1bSJed Brown         /* de = be * ae; */
554c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
555c4762a1bSJed Brown         ge = coef*be*(te - t0);
556c4762a1bSJed Brown 
557c4762a1bSJed Brown         ts = x[j-1][i];
558c4762a1bSJed Brown         as = 0.5*(t0 + ts);
559c4762a1bSJed Brown         bs = PetscPowScalar(as,bm1);
560c4762a1bSJed Brown         /* ds = bs * as; */
561c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
562c4762a1bSJed Brown         gs = coef*bs*(t0 - ts);
563c4762a1bSJed Brown 
564c4762a1bSJed Brown         v[0] = -hxdhy*(ds - gs);                             col[0].j = j-1;       col[0].i = i;
565c4762a1bSJed Brown         v[1] = -hydhx*(dw - gw);                             col[1].j = j;         col[1].i = i-1;
566c4762a1bSJed Brown         v[2] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge);  col[2].j = row.j = j; col[2].i = row.i = i;
567c4762a1bSJed Brown         v[3] = -hydhx*(de + ge);                             col[3].j = j;         col[3].i = i+1;
5689566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
569c4762a1bSJed Brown 
570c4762a1bSJed Brown       }
571c4762a1bSJed Brown     }
572c4762a1bSJed Brown   }
5739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
5749566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,localX,&x));
5759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
5769566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localX));
577c4762a1bSJed Brown   if (jac != B) {
5789566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
5799566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
580c4762a1bSJed Brown   }
581c4762a1bSJed Brown 
5829566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops((41.0 + 8.0*POWFLOP)*xm*ym));
583c4762a1bSJed Brown   PetscFunctionReturn(0);
584c4762a1bSJed Brown }
585c4762a1bSJed Brown 
586c4762a1bSJed Brown /*TEST
587c4762a1bSJed Brown 
588c4762a1bSJed Brown    test:
58941ba4c6cSHeeho Park       suffix: 1
590c4762a1bSJed Brown       args: -pc_type mg -ksp_type fgmres -da_refine 2 -pc_mg_galerkin pmat -snes_view
591c4762a1bSJed Brown       requires: !single
592c4762a1bSJed Brown 
59341ba4c6cSHeeho Park    test:
59441ba4c6cSHeeho Park       suffix: 2
59541ba4c6cSHeeho Park       args: -pc_type mg -ksp_type fgmres -da_refine 2 -pc_mg_galerkin pmat -snes_view -snes_type newtontrdc
59641ba4c6cSHeeho Park       requires: !single
59741ba4c6cSHeeho Park 
598c4762a1bSJed Brown TEST*/
599