xref: /petsc/src/snes/tutorials/ex18.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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 
609566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* set problem parameters */
63c4762a1bSJed Brown   user.tleft  = 1.0;
64c4762a1bSJed Brown   user.tright = 0.1;
65c4762a1bSJed Brown   user.beta   = 2.5;
669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-tleft",&user.tleft,NULL));
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-tright",&user.tright,NULL));
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-beta",&user.beta,NULL));
69c4762a1bSJed Brown   user.bm1    = user.beta - 1.0;
70c4762a1bSJed Brown   user.coef   = user.beta/2.0;
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /*
73c4762a1bSJed Brown       Create the multilevel DM data structure
74c4762a1bSJed Brown   */
759566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   /*
78c4762a1bSJed Brown       Set the DMDA (grid structure) for the grids.
79c4762a1bSJed Brown   */
809566063dSJacob 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));
819566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
829566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
839566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da,&user));
849566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes,(DM)da));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   /*
87c4762a1bSJed Brown      Create the nonlinear solver, and tell it the functions to use
88c4762a1bSJed Brown   */
899566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes,NULL,FormFunction,&user));
909566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes,NULL,NULL,FormJacobian,&user));
919566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
929566063dSJacob Faibussowitsch   PetscCall(SNESSetComputeInitialGuess(snes,FormInitialGuess,NULL));
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,NULL,NULL));
959566063dSJacob Faibussowitsch   PetscCall(SNESGetIterationNumber(snes,&its));
969566063dSJacob Faibussowitsch   PetscCall(SNESGetLinearSolveIterations(snes,&lits));
97c4762a1bSJed Brown   litspit = ((PetscReal)lits)/((PetscReal)its);
98*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of SNES iterations = %" PetscInt_FMT "\n",its));
99*63a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Number of Linear iterations = %" PetscInt_FMT "\n",lits));
1009566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Average Linear its / SNES = %e\n",(double)litspit));
101c4762a1bSJed Brown 
1029566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1039566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1049566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
105b122ec5aSJacob Faibussowitsch   return 0;
106c4762a1bSJed Brown }
107c4762a1bSJed Brown /* --------------------  Form initial approximation ----------------- */
108c4762a1bSJed Brown PetscErrorCode FormInitialGuess(SNES snes,Vec X,void *ctx)
109c4762a1bSJed Brown {
110c4762a1bSJed Brown   AppCtx         *user;
111c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym;
112c4762a1bSJed Brown   PetscReal      tleft;
113c4762a1bSJed Brown   PetscScalar    **x;
114c4762a1bSJed Brown   DM             da;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   PetscFunctionBeginUser;
1179566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&da));
1189566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(da,&user));
119c4762a1bSJed Brown   tleft = user->tleft;
120c4762a1bSJed Brown   /* Get ghost points */
1219566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0));
1229566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,X,&x));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   /* Compute initial guess */
125c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
126c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
127c4762a1bSJed Brown       x[j][i] = tleft;
128c4762a1bSJed Brown     }
129c4762a1bSJed Brown   }
1309566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,X,&x));
131c4762a1bSJed Brown   PetscFunctionReturn(0);
132c4762a1bSJed Brown }
133c4762a1bSJed Brown /* --------------------  Evaluate Function F(x) --------------------- */
134c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ptr)
135c4762a1bSJed Brown {
136c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
137c4762a1bSJed Brown   PetscInt       i,j,mx,my,xs,ys,xm,ym;
138c4762a1bSJed Brown   PetscScalar    zero = 0.0,one = 1.0;
139c4762a1bSJed Brown   PetscScalar    hx,hy,hxdhy,hydhx;
140c4762a1bSJed 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;
141c4762a1bSJed Brown   PetscScalar    tleft,tright,beta;
142c4762a1bSJed Brown   PetscScalar    **x,**f;
143c4762a1bSJed Brown   Vec            localX;
144c4762a1bSJed Brown   DM             da;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   PetscFunctionBeginUser;
1479566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&da));
1489566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localX));
1499566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,NULL,&mx,&my,0,0,0,0,0,0,0,0,0,0));
150c4762a1bSJed Brown   hx    = one/(PetscReal)(mx-1);  hy    = one/(PetscReal)(my-1);
151c4762a1bSJed Brown   hxdhy = hx/hy;               hydhx = hy/hx;
152c4762a1bSJed Brown   tleft = user->tleft;         tright = user->tright;
153c4762a1bSJed Brown   beta  = user->beta;
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* Get ghost points */
1569566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
1579566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
1589566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0));
1599566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,localX,&x));
1609566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,F,&f));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /* Evaluate function */
163c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
164c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
165c4762a1bSJed Brown       t0 = x[j][i];
166c4762a1bSJed Brown 
167c4762a1bSJed Brown       if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
168c4762a1bSJed Brown 
169c4762a1bSJed Brown         /* general interior volume */
170c4762a1bSJed Brown 
171c4762a1bSJed Brown         tw = x[j][i-1];
172c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
173c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
174c4762a1bSJed Brown         fw = dw*(t0 - tw);
175c4762a1bSJed Brown 
176c4762a1bSJed Brown         te = x[j][i+1];
177c4762a1bSJed Brown         ae = 0.5*(t0 + te);
178c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
179c4762a1bSJed Brown         fe = de*(te - t0);
180c4762a1bSJed Brown 
181c4762a1bSJed Brown         ts = x[j-1][i];
182c4762a1bSJed Brown         as = 0.5*(t0 + ts);
183c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
184c4762a1bSJed Brown         fs = ds*(t0 - ts);
185c4762a1bSJed Brown 
186c4762a1bSJed Brown         tn = x[j+1][i];
187c4762a1bSJed Brown         an = 0.5*(t0 + tn);
188c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
189c4762a1bSJed Brown         fn = dn*(tn - t0);
190c4762a1bSJed Brown 
191c4762a1bSJed Brown       } else if (i == 0) {
192c4762a1bSJed Brown 
193c4762a1bSJed Brown         /* left-hand boundary */
194c4762a1bSJed Brown         tw = tleft;
195c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
196c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
197c4762a1bSJed Brown         fw = dw*(t0 - tw);
198c4762a1bSJed Brown 
199c4762a1bSJed Brown         te = x[j][i+1];
200c4762a1bSJed Brown         ae = 0.5*(t0 + te);
201c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
202c4762a1bSJed Brown         fe = de*(te - t0);
203c4762a1bSJed Brown 
204c4762a1bSJed Brown         if (j > 0) {
205c4762a1bSJed Brown           ts = x[j-1][i];
206c4762a1bSJed Brown           as = 0.5*(t0 + ts);
207c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
208c4762a1bSJed Brown           fs = ds*(t0 - ts);
209c4762a1bSJed Brown         } else fs = zero;
210c4762a1bSJed Brown 
211c4762a1bSJed Brown         if (j < my-1) {
212c4762a1bSJed Brown           tn = x[j+1][i];
213c4762a1bSJed Brown           an = 0.5*(t0 + tn);
214c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
215c4762a1bSJed Brown           fn = dn*(tn - t0);
216c4762a1bSJed Brown         } else fn = zero;
217c4762a1bSJed Brown 
218c4762a1bSJed Brown       } else if (i == mx-1) {
219c4762a1bSJed Brown 
220c4762a1bSJed Brown         /* right-hand boundary */
221c4762a1bSJed Brown         tw = x[j][i-1];
222c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
223c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
224c4762a1bSJed Brown         fw = dw*(t0 - tw);
225c4762a1bSJed Brown 
226c4762a1bSJed Brown         te = tright;
227c4762a1bSJed Brown         ae = 0.5*(t0 + te);
228c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
229c4762a1bSJed Brown         fe = de*(te - t0);
230c4762a1bSJed Brown 
231c4762a1bSJed Brown         if (j > 0) {
232c4762a1bSJed Brown           ts = x[j-1][i];
233c4762a1bSJed Brown           as = 0.5*(t0 + ts);
234c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
235c4762a1bSJed Brown           fs = ds*(t0 - ts);
236c4762a1bSJed Brown         } else fs = zero;
237c4762a1bSJed Brown 
238c4762a1bSJed Brown         if (j < my-1) {
239c4762a1bSJed Brown           tn = x[j+1][i];
240c4762a1bSJed Brown           an = 0.5*(t0 + tn);
241c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
242c4762a1bSJed Brown           fn = dn*(tn - t0);
243c4762a1bSJed Brown         } else fn = zero;
244c4762a1bSJed Brown 
245c4762a1bSJed Brown       } else if (j == 0) {
246c4762a1bSJed Brown 
247c4762a1bSJed Brown         /* bottom boundary,and i <> 0 or mx-1 */
248c4762a1bSJed Brown         tw = x[j][i-1];
249c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
250c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
251c4762a1bSJed Brown         fw = dw*(t0 - tw);
252c4762a1bSJed Brown 
253c4762a1bSJed Brown         te = x[j][i+1];
254c4762a1bSJed Brown         ae = 0.5*(t0 + te);
255c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
256c4762a1bSJed Brown         fe = de*(te - t0);
257c4762a1bSJed Brown 
258c4762a1bSJed Brown         fs = zero;
259c4762a1bSJed Brown 
260c4762a1bSJed Brown         tn = x[j+1][i];
261c4762a1bSJed Brown         an = 0.5*(t0 + tn);
262c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
263c4762a1bSJed Brown         fn = dn*(tn - t0);
264c4762a1bSJed Brown 
265c4762a1bSJed Brown       } else if (j == my-1) {
266c4762a1bSJed Brown 
267c4762a1bSJed Brown         /* top boundary,and i <> 0 or mx-1 */
268c4762a1bSJed Brown         tw = x[j][i-1];
269c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
270c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
271c4762a1bSJed Brown         fw = dw*(t0 - tw);
272c4762a1bSJed Brown 
273c4762a1bSJed Brown         te = x[j][i+1];
274c4762a1bSJed Brown         ae = 0.5*(t0 + te);
275c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
276c4762a1bSJed Brown         fe = de*(te - t0);
277c4762a1bSJed Brown 
278c4762a1bSJed Brown         ts = x[j-1][i];
279c4762a1bSJed Brown         as = 0.5*(t0 + ts);
280c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
281c4762a1bSJed Brown         fs = ds*(t0 - ts);
282c4762a1bSJed Brown 
283c4762a1bSJed Brown         fn = zero;
284c4762a1bSJed Brown 
285c4762a1bSJed Brown       }
286c4762a1bSJed Brown 
287c4762a1bSJed Brown       f[j][i] = -hydhx*(fe-fw) - hxdhy*(fn-fs);
288c4762a1bSJed Brown 
289c4762a1bSJed Brown     }
290c4762a1bSJed Brown   }
2919566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,localX,&x));
2929566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,F,&f));
2939566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localX));
2949566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops((22.0 + 4.0*POWFLOP)*ym*xm));
295c4762a1bSJed Brown   PetscFunctionReturn(0);
296c4762a1bSJed Brown }
297c4762a1bSJed Brown /* --------------------  Evaluate Jacobian F(x) --------------------- */
298c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X,Mat jac,Mat B,void *ptr)
299c4762a1bSJed Brown {
300c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
301c4762a1bSJed Brown   PetscInt       i,j,mx,my,xs,ys,xm,ym;
302c4762a1bSJed Brown   PetscScalar    one = 1.0,hx,hy,hxdhy,hydhx,t0,tn,ts,te,tw;
303c4762a1bSJed Brown   PetscScalar    dn,ds,de,dw,an,as,ae,aw,bn,bs,be,bw,gn,gs,ge,gw;
304c4762a1bSJed Brown   PetscScalar    tleft,tright,beta,bm1,coef;
305c4762a1bSJed Brown   PetscScalar    v[5],**x;
306c4762a1bSJed Brown   Vec            localX;
307c4762a1bSJed Brown   MatStencil     col[5],row;
308c4762a1bSJed Brown   DM             da;
309c4762a1bSJed Brown 
310c4762a1bSJed Brown   PetscFunctionBeginUser;
3119566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&da));
3129566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localX));
3139566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,NULL,&mx,&my,0,0,0,0,0,0,0,0,0,0));
314c4762a1bSJed Brown   hx    = one/(PetscReal)(mx-1);  hy     = one/(PetscReal)(my-1);
315c4762a1bSJed Brown   hxdhy = hx/hy;               hydhx  = hy/hx;
316c4762a1bSJed Brown   tleft = user->tleft;         tright = user->tright;
317c4762a1bSJed Brown   beta  = user->beta;          bm1    = user->bm1;          coef = user->coef;
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   /* Get ghost points */
3209566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
3219566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
3229566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,0,&xm,&ym,0));
3239566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,localX,&x));
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   /* Evaluate Jacobian of function */
326c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
327c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
328c4762a1bSJed Brown       t0 = x[j][i];
329c4762a1bSJed Brown 
330c4762a1bSJed Brown       if (i > 0 && i < mx-1 && j > 0 && j < my-1) {
331c4762a1bSJed Brown 
332c4762a1bSJed Brown         /* general interior volume */
333c4762a1bSJed Brown 
334c4762a1bSJed Brown         tw = x[j][i-1];
335c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
336c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
337c4762a1bSJed Brown         /* dw = bw * aw */
338c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
339c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
340c4762a1bSJed Brown 
341c4762a1bSJed Brown         te = x[j][i+1];
342c4762a1bSJed Brown         ae = 0.5*(t0 + te);
343c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
344c4762a1bSJed Brown         /* de = be * ae; */
345c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
346c4762a1bSJed Brown         ge = coef*be*(te - t0);
347c4762a1bSJed Brown 
348c4762a1bSJed Brown         ts = x[j-1][i];
349c4762a1bSJed Brown         as = 0.5*(t0 + ts);
350c4762a1bSJed Brown         bs = PetscPowScalar(as,bm1);
351c4762a1bSJed Brown         /* ds = bs * as; */
352c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
353c4762a1bSJed Brown         gs = coef*bs*(t0 - ts);
354c4762a1bSJed Brown 
355c4762a1bSJed Brown         tn = x[j+1][i];
356c4762a1bSJed Brown         an = 0.5*(t0 + tn);
357c4762a1bSJed Brown         bn = PetscPowScalar(an,bm1);
358c4762a1bSJed Brown         /* dn = bn * an; */
359c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
360c4762a1bSJed Brown         gn = coef*bn*(tn - t0);
361c4762a1bSJed Brown 
362c4762a1bSJed Brown         v[0] = -hxdhy*(ds - gs);                                       col[0].j = j-1;       col[0].i = i;
363c4762a1bSJed Brown         v[1] = -hydhx*(dw - gw);                                       col[1].j = j;         col[1].i = i-1;
364c4762a1bSJed 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;
365c4762a1bSJed Brown         v[3] = -hydhx*(de + ge);                                       col[3].j = j;         col[3].i = i+1;
366c4762a1bSJed Brown         v[4] = -hxdhy*(dn + gn);                                       col[4].j = j+1;       col[4].i = i;
3679566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B,1,&row,5,col,v,INSERT_VALUES));
368c4762a1bSJed Brown 
369c4762a1bSJed Brown       } else if (i == 0) {
370c4762a1bSJed Brown 
371c4762a1bSJed Brown         /* left-hand boundary */
372c4762a1bSJed Brown         tw = tleft;
373c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
374c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
375c4762a1bSJed Brown         /* dw = bw * aw */
376c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
377c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
378c4762a1bSJed Brown 
379c4762a1bSJed Brown         te = x[j][i + 1];
380c4762a1bSJed Brown         ae = 0.5*(t0 + te);
381c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
382c4762a1bSJed Brown         /* de = be * ae; */
383c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
384c4762a1bSJed Brown         ge = coef*be*(te - t0);
385c4762a1bSJed Brown 
386c4762a1bSJed Brown         /* left-hand bottom boundary */
387c4762a1bSJed Brown         if (j == 0) {
388c4762a1bSJed Brown 
389c4762a1bSJed Brown           tn = x[j+1][i];
390c4762a1bSJed Brown           an = 0.5*(t0 + tn);
391c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
392c4762a1bSJed Brown           /* dn = bn * an; */
393c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
394c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
395c4762a1bSJed Brown 
396c4762a1bSJed Brown           v[0] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[0].j = row.j = j; col[0].i = row.i = i;
397c4762a1bSJed Brown           v[1] = -hydhx*(de + ge);                            col[1].j = j;         col[1].i = i+1;
398c4762a1bSJed Brown           v[2] = -hxdhy*(dn + gn);                            col[2].j = j+1;       col[2].i = i;
3999566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
400c4762a1bSJed Brown 
401c4762a1bSJed Brown           /* left-hand interior boundary */
402c4762a1bSJed Brown         } else if (j < my-1) {
403c4762a1bSJed Brown 
404c4762a1bSJed Brown           ts = x[j-1][i];
405c4762a1bSJed Brown           as = 0.5*(t0 + ts);
406c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
407c4762a1bSJed Brown           /* ds = bs * as; */
408c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
409c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
410c4762a1bSJed Brown 
411c4762a1bSJed Brown           tn = x[j+1][i];
412c4762a1bSJed Brown           an = 0.5*(t0 + tn);
413c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
414c4762a1bSJed Brown           /* dn = bn * an; */
415c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
416c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
417c4762a1bSJed Brown 
418c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                                       col[0].j = j-1;       col[0].i = i;
419c4762a1bSJed 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;
420c4762a1bSJed Brown           v[2] = -hydhx*(de + ge);                                       col[2].j = j;         col[2].i = i+1;
421c4762a1bSJed Brown           v[3] = -hxdhy*(dn + gn);                                       col[3].j = j+1;       col[3].i = i;
4229566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
423c4762a1bSJed Brown           /* left-hand top boundary */
424c4762a1bSJed Brown         } else {
425c4762a1bSJed Brown 
426c4762a1bSJed Brown           ts = x[j-1][i];
427c4762a1bSJed Brown           as = 0.5*(t0 + ts);
428c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
429c4762a1bSJed Brown           /* ds = bs * as; */
430c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
431c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
432c4762a1bSJed Brown 
433c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                             col[0].j = j-1;       col[0].i = i;
434c4762a1bSJed Brown           v[1] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge);  col[1].j = row.j = j; col[1].i = row.i = i;
435c4762a1bSJed Brown           v[2] = -hydhx*(de + ge);                             col[2].j = j;         col[2].i = i+1;
4369566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
437c4762a1bSJed Brown         }
438c4762a1bSJed Brown 
439c4762a1bSJed Brown       } else if (i == mx-1) {
440c4762a1bSJed Brown 
441c4762a1bSJed Brown         /* right-hand boundary */
442c4762a1bSJed Brown         tw = x[j][i-1];
443c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
444c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
445c4762a1bSJed Brown         /* dw = bw * aw */
446c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
447c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
448c4762a1bSJed Brown 
449c4762a1bSJed Brown         te = tright;
450c4762a1bSJed Brown         ae = 0.5*(t0 + te);
451c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
452c4762a1bSJed Brown         /* de = be * ae; */
453c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
454c4762a1bSJed Brown         ge = coef*be*(te - t0);
455c4762a1bSJed Brown 
456c4762a1bSJed Brown         /* right-hand bottom boundary */
457c4762a1bSJed Brown         if (j == 0) {
458c4762a1bSJed Brown 
459c4762a1bSJed Brown           tn = x[j+1][i];
460c4762a1bSJed Brown           an = 0.5*(t0 + tn);
461c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
462c4762a1bSJed Brown           /* dn = bn * an; */
463c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
464c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
465c4762a1bSJed Brown 
466c4762a1bSJed Brown           v[0] = -hydhx*(dw - gw);                            col[0].j = j;         col[0].i = i-1;
467c4762a1bSJed Brown           v[1] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[1].j = row.j = j; col[1].i = row.i = i;
468c4762a1bSJed Brown           v[2] = -hxdhy*(dn + gn);                            col[2].j = j+1;       col[2].i = i;
4699566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
470c4762a1bSJed Brown 
471c4762a1bSJed Brown           /* right-hand interior boundary */
472c4762a1bSJed Brown         } else if (j < my-1) {
473c4762a1bSJed Brown 
474c4762a1bSJed Brown           ts = x[j-1][i];
475c4762a1bSJed Brown           as = 0.5*(t0 + ts);
476c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
477c4762a1bSJed Brown           /* ds = bs * as; */
478c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
479c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
480c4762a1bSJed Brown 
481c4762a1bSJed Brown           tn = x[j+1][i];
482c4762a1bSJed Brown           an = 0.5*(t0 + tn);
483c4762a1bSJed Brown           bn = PetscPowScalar(an,bm1);
484c4762a1bSJed Brown           /* dn = bn * an; */
485c4762a1bSJed Brown           dn = PetscPowScalar(an,beta);
486c4762a1bSJed Brown           gn = coef*bn*(tn - t0);
487c4762a1bSJed Brown 
488c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                                       col[0].j = j-1;       col[0].i = i;
489c4762a1bSJed Brown           v[1] = -hydhx*(dw - gw);                                       col[1].j = j;         col[1].i = i-1;
490c4762a1bSJed 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;
491c4762a1bSJed Brown           v[3] = -hxdhy*(dn + gn);                                       col[3].j = j+1;       col[3].i = i;
4929566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
493c4762a1bSJed Brown         /* right-hand top boundary */
494c4762a1bSJed Brown         } else {
495c4762a1bSJed Brown 
496c4762a1bSJed Brown           ts = x[j-1][i];
497c4762a1bSJed Brown           as = 0.5*(t0 + ts);
498c4762a1bSJed Brown           bs = PetscPowScalar(as,bm1);
499c4762a1bSJed Brown           /* ds = bs * as; */
500c4762a1bSJed Brown           ds = PetscPowScalar(as,beta);
501c4762a1bSJed Brown           gs = coef*bs*(t0 - ts);
502c4762a1bSJed Brown 
503c4762a1bSJed Brown           v[0] = -hxdhy*(ds - gs);                             col[0].j = j-1;       col[0].i = i;
504c4762a1bSJed Brown           v[1] = -hydhx*(dw - gw);                             col[1].j = j;         col[1].i = i-1;
505c4762a1bSJed Brown           v[2] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge);  col[2].j = row.j = j; col[2].i = row.i = i;
5069566063dSJacob Faibussowitsch           PetscCall(MatSetValuesStencil(B,1,&row,3,col,v,INSERT_VALUES));
507c4762a1bSJed Brown         }
508c4762a1bSJed Brown 
509c4762a1bSJed Brown         /* bottom boundary,and i <> 0 or mx-1 */
510c4762a1bSJed Brown       } else if (j == 0) {
511c4762a1bSJed Brown 
512c4762a1bSJed Brown         tw = x[j][i-1];
513c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
514c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
515c4762a1bSJed Brown         /* dw = bw * aw */
516c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
517c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
518c4762a1bSJed Brown 
519c4762a1bSJed Brown         te = x[j][i+1];
520c4762a1bSJed Brown         ae = 0.5*(t0 + te);
521c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
522c4762a1bSJed Brown         /* de = be * ae; */
523c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
524c4762a1bSJed Brown         ge = coef*be*(te - t0);
525c4762a1bSJed Brown 
526c4762a1bSJed Brown         tn = x[j+1][i];
527c4762a1bSJed Brown         an = 0.5*(t0 + tn);
528c4762a1bSJed Brown         bn = PetscPowScalar(an,bm1);
529c4762a1bSJed Brown         /* dn = bn * an; */
530c4762a1bSJed Brown         dn = PetscPowScalar(an,beta);
531c4762a1bSJed Brown         gn = coef*bn*(tn - t0);
532c4762a1bSJed Brown 
533c4762a1bSJed Brown         v[0] = -hydhx*(dw - gw);                            col[0].j = j;         col[0].i = i-1;
534c4762a1bSJed Brown         v[1] = hxdhy*(dn - gn) + hydhx*(dw + de + gw - ge); col[1].j = row.j = j; col[1].i = row.i = i;
535c4762a1bSJed Brown         v[2] = -hydhx*(de + ge);                            col[2].j = j;         col[2].i = i+1;
536c4762a1bSJed Brown         v[3] = -hxdhy*(dn + gn);                            col[3].j = j+1;       col[3].i = i;
5379566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
538c4762a1bSJed Brown 
539c4762a1bSJed Brown         /* top boundary,and i <> 0 or mx-1 */
540c4762a1bSJed Brown       } else if (j == my-1) {
541c4762a1bSJed Brown 
542c4762a1bSJed Brown         tw = x[j][i-1];
543c4762a1bSJed Brown         aw = 0.5*(t0 + tw);
544c4762a1bSJed Brown         bw = PetscPowScalar(aw,bm1);
545c4762a1bSJed Brown         /* dw = bw * aw */
546c4762a1bSJed Brown         dw = PetscPowScalar(aw,beta);
547c4762a1bSJed Brown         gw = coef*bw*(t0 - tw);
548c4762a1bSJed Brown 
549c4762a1bSJed Brown         te = x[j][i+1];
550c4762a1bSJed Brown         ae = 0.5*(t0 + te);
551c4762a1bSJed Brown         be = PetscPowScalar(ae,bm1);
552c4762a1bSJed Brown         /* de = be * ae; */
553c4762a1bSJed Brown         de = PetscPowScalar(ae,beta);
554c4762a1bSJed Brown         ge = coef*be*(te - t0);
555c4762a1bSJed Brown 
556c4762a1bSJed Brown         ts = x[j-1][i];
557c4762a1bSJed Brown         as = 0.5*(t0 + ts);
558c4762a1bSJed Brown         bs = PetscPowScalar(as,bm1);
559c4762a1bSJed Brown         /* ds = bs * as; */
560c4762a1bSJed Brown         ds = PetscPowScalar(as,beta);
561c4762a1bSJed Brown         gs = coef*bs*(t0 - ts);
562c4762a1bSJed Brown 
563c4762a1bSJed Brown         v[0] = -hxdhy*(ds - gs);                             col[0].j = j-1;       col[0].i = i;
564c4762a1bSJed Brown         v[1] = -hydhx*(dw - gw);                             col[1].j = j;         col[1].i = i-1;
565c4762a1bSJed Brown         v[2] = hxdhy*(ds + gs) + hydhx*(dw + de + gw - ge);  col[2].j = row.j = j; col[2].i = row.i = i;
566c4762a1bSJed Brown         v[3] = -hydhx*(de + ge);                             col[3].j = j;         col[3].i = i+1;
5679566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B,1,&row,4,col,v,INSERT_VALUES));
568c4762a1bSJed Brown 
569c4762a1bSJed Brown       }
570c4762a1bSJed Brown     }
571c4762a1bSJed Brown   }
5729566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
5739566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,localX,&x));
5749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
5759566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localX));
576c4762a1bSJed Brown   if (jac != B) {
5779566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
5789566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
579c4762a1bSJed Brown   }
580c4762a1bSJed Brown 
5819566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops((41.0 + 8.0*POWFLOP)*xm*ym));
582c4762a1bSJed Brown   PetscFunctionReturn(0);
583c4762a1bSJed Brown }
584c4762a1bSJed Brown 
585c4762a1bSJed Brown /*TEST
586c4762a1bSJed Brown 
587c4762a1bSJed Brown    test:
58841ba4c6cSHeeho Park       suffix: 1
589c4762a1bSJed Brown       args: -pc_type mg -ksp_type fgmres -da_refine 2 -pc_mg_galerkin pmat -snes_view
590c4762a1bSJed Brown       requires: !single
591c4762a1bSJed Brown 
59241ba4c6cSHeeho Park    test:
59341ba4c6cSHeeho Park       suffix: 2
59441ba4c6cSHeeho Park       args: -pc_type mg -ksp_type fgmres -da_refine 2 -pc_mg_galerkin pmat -snes_view -snes_type newtontrdc
59541ba4c6cSHeeho Park       requires: !single
59641ba4c6cSHeeho Park 
597c4762a1bSJed Brown TEST*/
598