xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/reaction_diffusion.c (revision 60f0b76ed2293387b137bb632a681b21c2f7fb7d)
1*60f0b76eSHong Zhang #include "reaction_diffusion.h"
2*60f0b76eSHong Zhang #include <petscdm.h>
3*60f0b76eSHong Zhang #include <petscdmda.h>
4*60f0b76eSHong Zhang 
5*60f0b76eSHong Zhang /*F
6*60f0b76eSHong Zhang      This example is taken from the book, Numerical Solution of Time-Dependent Advection-Diffusion-Reaction Equations by
7*60f0b76eSHong Zhang       W. Hundsdorf and J.G. Verwer,  Page 21, Pattern Formation with Reaction-Diffusion Equations
8*60f0b76eSHong Zhang \begin{eqnarray*}
9*60f0b76eSHong Zhang         u_t = D_1 (u_{xx} + u_{yy})  - u*v^2 + \gamma(1 -u)           \\
10*60f0b76eSHong Zhang         v_t = D_2 (v_{xx} + v_{yy})  + u*v^2 - (\gamma + \kappa)v
11*60f0b76eSHong Zhang \end{eqnarray*}
12*60f0b76eSHong Zhang     Unlike in the book this uses periodic boundary conditions instead of Neumann
13*60f0b76eSHong Zhang     (since they are easier for finite differences).
14*60f0b76eSHong Zhang F*/
15*60f0b76eSHong Zhang 
16*60f0b76eSHong Zhang /*
17*60f0b76eSHong Zhang    RHSFunction - Evaluates nonlinear function, F(x).
18*60f0b76eSHong Zhang 
19*60f0b76eSHong Zhang    Input Parameters:
20*60f0b76eSHong Zhang .  ts - the TS context
21*60f0b76eSHong Zhang .  X - input vector
22*60f0b76eSHong Zhang .  ptr - optional user-defined context, as set by TSSetRHSFunction()
23*60f0b76eSHong Zhang 
24*60f0b76eSHong Zhang    Output Parameter:
25*60f0b76eSHong Zhang .  F - function vector
26*60f0b76eSHong Zhang  */
27*60f0b76eSHong Zhang PetscErrorCode RHSFunction(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr)
28*60f0b76eSHong Zhang {
29*60f0b76eSHong Zhang   AppCtx         *appctx = (AppCtx*)ptr;
30*60f0b76eSHong Zhang   DM             da;
31*60f0b76eSHong Zhang   PetscErrorCode ierr;
32*60f0b76eSHong Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
33*60f0b76eSHong Zhang   PetscReal      hx,hy,sx,sy;
34*60f0b76eSHong Zhang   PetscScalar    uc,uxx,uyy,vc,vxx,vyy;
35*60f0b76eSHong Zhang   Field          **u,**f;
36*60f0b76eSHong Zhang   Vec            localU;
37*60f0b76eSHong Zhang 
38*60f0b76eSHong Zhang   PetscFunctionBegin;
39*60f0b76eSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
40*60f0b76eSHong Zhang   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
41*60f0b76eSHong Zhang   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
42*60f0b76eSHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
43*60f0b76eSHong Zhang   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
44*60f0b76eSHong Zhang 
45*60f0b76eSHong Zhang   /*
46*60f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
47*60f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
48*60f0b76eSHong Zhang      By placing code between these two statements, computations can be
49*60f0b76eSHong Zhang      done while messages are in transition.
50*60f0b76eSHong Zhang   */
51*60f0b76eSHong Zhang   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
52*60f0b76eSHong Zhang   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
53*60f0b76eSHong Zhang 
54*60f0b76eSHong Zhang   /*
55*60f0b76eSHong Zhang      Get pointers to vector data
56*60f0b76eSHong Zhang   */
57*60f0b76eSHong Zhang   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
58*60f0b76eSHong Zhang   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
59*60f0b76eSHong Zhang 
60*60f0b76eSHong Zhang   /*
61*60f0b76eSHong Zhang      Get local grid boundaries
62*60f0b76eSHong Zhang   */
63*60f0b76eSHong Zhang   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
64*60f0b76eSHong Zhang 
65*60f0b76eSHong Zhang   /*
66*60f0b76eSHong Zhang      Compute function over the locally owned part of the grid
67*60f0b76eSHong Zhang   */
68*60f0b76eSHong Zhang   for (j=ys; j<ys+ym; j++) {
69*60f0b76eSHong Zhang     for (i=xs; i<xs+xm; i++) {
70*60f0b76eSHong Zhang       uc        = u[j][i].u;
71*60f0b76eSHong Zhang       uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
72*60f0b76eSHong Zhang       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
73*60f0b76eSHong Zhang       vc        = u[j][i].v;
74*60f0b76eSHong Zhang       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
75*60f0b76eSHong Zhang       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
76*60f0b76eSHong Zhang       f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);
77*60f0b76eSHong Zhang       f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc;
78*60f0b76eSHong Zhang     }
79*60f0b76eSHong Zhang   }
80*60f0b76eSHong Zhang   ierr = PetscLogFlops(16.0*xm*ym);CHKERRQ(ierr);
81*60f0b76eSHong Zhang 
82*60f0b76eSHong Zhang   /*
83*60f0b76eSHong Zhang      Restore vectors
84*60f0b76eSHong Zhang   */
85*60f0b76eSHong Zhang   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
86*60f0b76eSHong Zhang   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
87*60f0b76eSHong Zhang   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
88*60f0b76eSHong Zhang   PetscFunctionReturn(0);
89*60f0b76eSHong Zhang }
90*60f0b76eSHong Zhang 
91*60f0b76eSHong Zhang PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat BB,void *ctx)
92*60f0b76eSHong Zhang {
93*60f0b76eSHong Zhang   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
94*60f0b76eSHong Zhang   DM             da;
95*60f0b76eSHong Zhang   PetscErrorCode ierr;
96*60f0b76eSHong Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
97*60f0b76eSHong Zhang   PetscReal      hx,hy,sx,sy;
98*60f0b76eSHong Zhang   PetscScalar    uc,vc;
99*60f0b76eSHong Zhang   Field          **u;
100*60f0b76eSHong Zhang   Vec            localU;
101*60f0b76eSHong Zhang   MatStencil     stencil[6],rowstencil;
102*60f0b76eSHong Zhang   PetscScalar    entries[6];
103*60f0b76eSHong Zhang 
104*60f0b76eSHong Zhang   PetscFunctionBegin;
105*60f0b76eSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
106*60f0b76eSHong Zhang   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
107*60f0b76eSHong Zhang   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
108*60f0b76eSHong Zhang 
109*60f0b76eSHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
110*60f0b76eSHong Zhang   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
111*60f0b76eSHong Zhang 
112*60f0b76eSHong Zhang   /*
113*60f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
114*60f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
115*60f0b76eSHong Zhang      By placing code between these two statements, computations can be
116*60f0b76eSHong Zhang      done while messages are in transition.
117*60f0b76eSHong Zhang   */
118*60f0b76eSHong Zhang   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
119*60f0b76eSHong Zhang   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
120*60f0b76eSHong Zhang 
121*60f0b76eSHong Zhang   /*
122*60f0b76eSHong Zhang      Get pointers to vector data
123*60f0b76eSHong Zhang   */
124*60f0b76eSHong Zhang   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
125*60f0b76eSHong Zhang 
126*60f0b76eSHong Zhang   /*
127*60f0b76eSHong Zhang      Get local grid boundaries
128*60f0b76eSHong Zhang   */
129*60f0b76eSHong Zhang   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
130*60f0b76eSHong Zhang 
131*60f0b76eSHong Zhang   stencil[0].k = 0;
132*60f0b76eSHong Zhang   stencil[1].k = 0;
133*60f0b76eSHong Zhang   stencil[2].k = 0;
134*60f0b76eSHong Zhang   stencil[3].k = 0;
135*60f0b76eSHong Zhang   stencil[4].k = 0;
136*60f0b76eSHong Zhang   stencil[5].k = 0;
137*60f0b76eSHong Zhang   rowstencil.k = 0;
138*60f0b76eSHong Zhang   /*
139*60f0b76eSHong Zhang      Compute function over the locally owned part of the grid
140*60f0b76eSHong Zhang   */
141*60f0b76eSHong Zhang   for (j=ys; j<ys+ym; j++) {
142*60f0b76eSHong Zhang     stencil[0].j = j-1;
143*60f0b76eSHong Zhang     stencil[1].j = j+1;
144*60f0b76eSHong Zhang     stencil[2].j = j;
145*60f0b76eSHong Zhang     stencil[3].j = j;
146*60f0b76eSHong Zhang     stencil[4].j = j;
147*60f0b76eSHong Zhang     stencil[5].j = j;
148*60f0b76eSHong Zhang     rowstencil.k = 0; rowstencil.j = j;
149*60f0b76eSHong Zhang     for (i=xs; i<xs+xm; i++) {
150*60f0b76eSHong Zhang       uc = u[j][i].u;
151*60f0b76eSHong Zhang       vc = u[j][i].v;
152*60f0b76eSHong Zhang 
153*60f0b76eSHong Zhang       /*      uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
154*60f0b76eSHong Zhang       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
155*60f0b76eSHong Zhang 
156*60f0b76eSHong Zhang       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
157*60f0b76eSHong Zhang       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
158*60f0b76eSHong Zhang        f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
159*60f0b76eSHong Zhang 
160*60f0b76eSHong Zhang       stencil[0].i = i; stencil[0].c = 0; entries[0] = appctx->D1*sy;
161*60f0b76eSHong Zhang       stencil[1].i = i; stencil[1].c = 0; entries[1] = appctx->D1*sy;
162*60f0b76eSHong Zhang       stencil[2].i = i-1; stencil[2].c = 0; entries[2] = appctx->D1*sx;
163*60f0b76eSHong Zhang       stencil[3].i = i+1; stencil[3].c = 0; entries[3] = appctx->D1*sx;
164*60f0b76eSHong Zhang       stencil[4].i = i; stencil[4].c = 0; entries[4] = -2.0*appctx->D1*(sx + sy) - vc*vc - appctx->gamma;
165*60f0b76eSHong Zhang       stencil[5].i = i; stencil[5].c = 1; entries[5] = -2.0*uc*vc;
166*60f0b76eSHong Zhang       rowstencil.i = i; rowstencil.c = 0;
167*60f0b76eSHong Zhang 
168*60f0b76eSHong Zhang       ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
169*60f0b76eSHong Zhang       if (appctx->aijpc) {
170*60f0b76eSHong Zhang         ierr = MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
171*60f0b76eSHong Zhang       }
172*60f0b76eSHong Zhang       stencil[0].c = 1; entries[0] = appctx->D2*sy;
173*60f0b76eSHong Zhang       stencil[1].c = 1; entries[1] = appctx->D2*sy;
174*60f0b76eSHong Zhang       stencil[2].c = 1; entries[2] = appctx->D2*sx;
175*60f0b76eSHong Zhang       stencil[3].c = 1; entries[3] = appctx->D2*sx;
176*60f0b76eSHong Zhang       stencil[4].c = 1; entries[4] = -2.0*appctx->D2*(sx + sy) + 2.0*uc*vc - appctx->gamma - appctx->kappa;
177*60f0b76eSHong Zhang       stencil[5].c = 0; entries[5] = vc*vc;
178*60f0b76eSHong Zhang       rowstencil.c = 1;
179*60f0b76eSHong Zhang 
180*60f0b76eSHong Zhang       ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
181*60f0b76eSHong Zhang       if (appctx->aijpc) {
182*60f0b76eSHong Zhang         ierr = MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
183*60f0b76eSHong Zhang       }
184*60f0b76eSHong Zhang       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
185*60f0b76eSHong Zhang     }
186*60f0b76eSHong Zhang   }
187*60f0b76eSHong Zhang 
188*60f0b76eSHong Zhang   /*
189*60f0b76eSHong Zhang      Restore vectors
190*60f0b76eSHong Zhang   */
191*60f0b76eSHong Zhang   ierr = PetscLogFlops(19.0*xm*ym);CHKERRQ(ierr);
192*60f0b76eSHong Zhang   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
193*60f0b76eSHong Zhang   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
194*60f0b76eSHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
195*60f0b76eSHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
196*60f0b76eSHong Zhang   ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
197*60f0b76eSHong Zhang   if (appctx->aijpc) {
198*60f0b76eSHong Zhang     ierr = MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
199*60f0b76eSHong Zhang     ierr = MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
200*60f0b76eSHong Zhang     ierr = MatSetOption(BB,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
201*60f0b76eSHong Zhang   }
202*60f0b76eSHong Zhang   PetscFunctionReturn(0);
203*60f0b76eSHong Zhang }
204*60f0b76eSHong Zhang 
205*60f0b76eSHong Zhang /*
206*60f0b76eSHong Zhang    IFunction - Evaluates implicit nonlinear function, xdot - F(x).
207*60f0b76eSHong Zhang 
208*60f0b76eSHong Zhang    Input Parameters:
209*60f0b76eSHong Zhang .  ts - the TS context
210*60f0b76eSHong Zhang .  U - input vector
211*60f0b76eSHong Zhang .  Udot - input vector
212*60f0b76eSHong Zhang .  ptr - optional user-defined context, as set by TSSetRHSFunction()
213*60f0b76eSHong Zhang 
214*60f0b76eSHong Zhang    Output Parameter:
215*60f0b76eSHong Zhang .  F - function vector
216*60f0b76eSHong Zhang  */
217*60f0b76eSHong Zhang PetscErrorCode IFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
218*60f0b76eSHong Zhang {
219*60f0b76eSHong Zhang   AppCtx         *appctx = (AppCtx*)ptr;
220*60f0b76eSHong Zhang   DM             da;
221*60f0b76eSHong Zhang   PetscErrorCode ierr;
222*60f0b76eSHong Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
223*60f0b76eSHong Zhang   PetscReal      hx,hy,sx,sy;
224*60f0b76eSHong Zhang   PetscScalar    uc,uxx,uyy,vc,vxx,vyy;
225*60f0b76eSHong Zhang   Field          **u,**f,**udot;
226*60f0b76eSHong Zhang   Vec            localU;
227*60f0b76eSHong Zhang 
228*60f0b76eSHong Zhang   PetscFunctionBegin;
229*60f0b76eSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
230*60f0b76eSHong Zhang   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
231*60f0b76eSHong Zhang   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
232*60f0b76eSHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
233*60f0b76eSHong Zhang   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
234*60f0b76eSHong Zhang 
235*60f0b76eSHong Zhang   /*
236*60f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
237*60f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
238*60f0b76eSHong Zhang      By placing code between these two statements, computations can be
239*60f0b76eSHong Zhang      done while messages are in transition.
240*60f0b76eSHong Zhang   */
241*60f0b76eSHong Zhang   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
242*60f0b76eSHong Zhang   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
243*60f0b76eSHong Zhang 
244*60f0b76eSHong Zhang   /*
245*60f0b76eSHong Zhang      Get pointers to vector data
246*60f0b76eSHong Zhang   */
247*60f0b76eSHong Zhang   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
248*60f0b76eSHong Zhang   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
249*60f0b76eSHong Zhang   ierr = DMDAVecGetArrayRead(da,Udot,&udot);CHKERRQ(ierr);
250*60f0b76eSHong Zhang 
251*60f0b76eSHong Zhang   /*
252*60f0b76eSHong Zhang      Get local grid boundaries
253*60f0b76eSHong Zhang   */
254*60f0b76eSHong Zhang   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
255*60f0b76eSHong Zhang 
256*60f0b76eSHong Zhang   /*
257*60f0b76eSHong Zhang      Compute function over the locally owned part of the grid
258*60f0b76eSHong Zhang   */
259*60f0b76eSHong Zhang   for (j=ys; j<ys+ym; j++) {
260*60f0b76eSHong Zhang     for (i=xs; i<xs+xm; i++) {
261*60f0b76eSHong Zhang       uc        = u[j][i].u;
262*60f0b76eSHong Zhang       uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
263*60f0b76eSHong Zhang       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
264*60f0b76eSHong Zhang       vc        = u[j][i].v;
265*60f0b76eSHong Zhang       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
266*60f0b76eSHong Zhang       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
267*60f0b76eSHong Zhang       f[j][i].u = udot[j][i].u - ( appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc));
268*60f0b76eSHong Zhang       f[j][i].v = udot[j][i].v - ( appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc);
269*60f0b76eSHong Zhang     }
270*60f0b76eSHong Zhang   }
271*60f0b76eSHong Zhang   ierr = PetscLogFlops(16.0*xm*ym);CHKERRQ(ierr);
272*60f0b76eSHong Zhang 
273*60f0b76eSHong Zhang   /*
274*60f0b76eSHong Zhang      Restore vectors
275*60f0b76eSHong Zhang   */
276*60f0b76eSHong Zhang   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
277*60f0b76eSHong Zhang   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
278*60f0b76eSHong Zhang   ierr = DMDAVecRestoreArrayRead(da,Udot,&udot);CHKERRQ(ierr);
279*60f0b76eSHong Zhang   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
280*60f0b76eSHong Zhang   PetscFunctionReturn(0);
281*60f0b76eSHong Zhang }
282*60f0b76eSHong Zhang 
283*60f0b76eSHong Zhang PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat BB,void *ctx)
284*60f0b76eSHong Zhang {
285*60f0b76eSHong Zhang   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
286*60f0b76eSHong Zhang   DM             da;
287*60f0b76eSHong Zhang   PetscErrorCode ierr;
288*60f0b76eSHong Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
289*60f0b76eSHong Zhang   PetscReal      hx,hy,sx,sy;
290*60f0b76eSHong Zhang   PetscScalar    uc,vc;
291*60f0b76eSHong Zhang   Field          **u;
292*60f0b76eSHong Zhang   Vec            localU;
293*60f0b76eSHong Zhang   MatStencil     stencil[6],rowstencil;
294*60f0b76eSHong Zhang   PetscScalar    entries[6];
295*60f0b76eSHong Zhang 
296*60f0b76eSHong Zhang   PetscFunctionBegin;
297*60f0b76eSHong Zhang   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
298*60f0b76eSHong Zhang   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
299*60f0b76eSHong Zhang   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
300*60f0b76eSHong Zhang 
301*60f0b76eSHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
302*60f0b76eSHong Zhang   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
303*60f0b76eSHong Zhang 
304*60f0b76eSHong Zhang   /*
305*60f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
306*60f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
307*60f0b76eSHong Zhang      By placing code between these two statements, computations can be
308*60f0b76eSHong Zhang      done while messages are in transition.
309*60f0b76eSHong Zhang   */
310*60f0b76eSHong Zhang   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
311*60f0b76eSHong Zhang   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
312*60f0b76eSHong Zhang 
313*60f0b76eSHong Zhang   /*
314*60f0b76eSHong Zhang      Get pointers to vector data
315*60f0b76eSHong Zhang   */
316*60f0b76eSHong Zhang   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
317*60f0b76eSHong Zhang 
318*60f0b76eSHong Zhang   /*
319*60f0b76eSHong Zhang      Get local grid boundaries
320*60f0b76eSHong Zhang   */
321*60f0b76eSHong Zhang   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
322*60f0b76eSHong Zhang 
323*60f0b76eSHong Zhang   stencil[0].k = 0;
324*60f0b76eSHong Zhang   stencil[1].k = 0;
325*60f0b76eSHong Zhang   stencil[2].k = 0;
326*60f0b76eSHong Zhang   stencil[3].k = 0;
327*60f0b76eSHong Zhang   stencil[4].k = 0;
328*60f0b76eSHong Zhang   stencil[5].k = 0;
329*60f0b76eSHong Zhang   rowstencil.k = 0;
330*60f0b76eSHong Zhang   /*
331*60f0b76eSHong Zhang      Compute function over the locally owned part of the grid
332*60f0b76eSHong Zhang   */
333*60f0b76eSHong Zhang   for (j=ys; j<ys+ym; j++) {
334*60f0b76eSHong Zhang 
335*60f0b76eSHong Zhang     stencil[0].j = j-1;
336*60f0b76eSHong Zhang     stencil[1].j = j+1;
337*60f0b76eSHong Zhang     stencil[2].j = j;
338*60f0b76eSHong Zhang     stencil[3].j = j;
339*60f0b76eSHong Zhang     stencil[4].j = j;
340*60f0b76eSHong Zhang     stencil[5].j = j;
341*60f0b76eSHong Zhang     rowstencil.k = 0; rowstencil.j = j;
342*60f0b76eSHong Zhang     for (i=xs; i<xs+xm; i++) {
343*60f0b76eSHong Zhang       uc = u[j][i].u;
344*60f0b76eSHong Zhang       vc = u[j][i].v;
345*60f0b76eSHong Zhang 
346*60f0b76eSHong Zhang       /*      uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
347*60f0b76eSHong Zhang       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
348*60f0b76eSHong Zhang 
349*60f0b76eSHong Zhang       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
350*60f0b76eSHong Zhang       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
351*60f0b76eSHong Zhang        f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
352*60f0b76eSHong Zhang 
353*60f0b76eSHong Zhang       stencil[0].i = i; stencil[0].c = 0; entries[0] = -appctx->D1*sy;
354*60f0b76eSHong Zhang       stencil[1].i = i; stencil[1].c = 0; entries[1] = -appctx->D1*sy;
355*60f0b76eSHong Zhang       stencil[2].i = i-1; stencil[2].c = 0; entries[2] = -appctx->D1*sx;
356*60f0b76eSHong Zhang       stencil[3].i = i+1; stencil[3].c = 0; entries[3] = -appctx->D1*sx;
357*60f0b76eSHong Zhang       stencil[4].i = i; stencil[4].c = 0; entries[4] = 2.0*appctx->D1*(sx + sy) + vc*vc + appctx->gamma + a;
358*60f0b76eSHong Zhang       stencil[5].i = i; stencil[5].c = 1; entries[5] = 2.0*uc*vc;
359*60f0b76eSHong Zhang       rowstencil.i = i; rowstencil.c = 0;
360*60f0b76eSHong Zhang 
361*60f0b76eSHong Zhang       ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
362*60f0b76eSHong Zhang       if (appctx->aijpc) {
363*60f0b76eSHong Zhang         ierr = MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
364*60f0b76eSHong Zhang       }
365*60f0b76eSHong Zhang       stencil[0].c = 1; entries[0] = -appctx->D2*sy;
366*60f0b76eSHong Zhang       stencil[1].c = 1; entries[1] = -appctx->D2*sy;
367*60f0b76eSHong Zhang       stencil[2].c = 1; entries[2] = -appctx->D2*sx;
368*60f0b76eSHong Zhang       stencil[3].c = 1; entries[3] = -appctx->D2*sx;
369*60f0b76eSHong Zhang       stencil[4].c = 1; entries[4] = 2.0*appctx->D2*(sx + sy) - 2.0*uc*vc + appctx->gamma + appctx->kappa + a;
370*60f0b76eSHong Zhang       stencil[5].c = 0; entries[5] = -vc*vc;
371*60f0b76eSHong Zhang       rowstencil.c = 1;
372*60f0b76eSHong Zhang 
373*60f0b76eSHong Zhang       ierr = MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
374*60f0b76eSHong Zhang       if (appctx->aijpc) {
375*60f0b76eSHong Zhang         ierr = MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES);CHKERRQ(ierr);
376*60f0b76eSHong Zhang       }
377*60f0b76eSHong Zhang       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
378*60f0b76eSHong Zhang     }
379*60f0b76eSHong Zhang   }
380*60f0b76eSHong Zhang 
381*60f0b76eSHong Zhang   /*
382*60f0b76eSHong Zhang      Restore vectors
383*60f0b76eSHong Zhang   */
384*60f0b76eSHong Zhang   ierr = PetscLogFlops(19.0*xm*ym);CHKERRQ(ierr);
385*60f0b76eSHong Zhang   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
386*60f0b76eSHong Zhang   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
387*60f0b76eSHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
388*60f0b76eSHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
389*60f0b76eSHong Zhang   ierr = MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
390*60f0b76eSHong Zhang   if (appctx->aijpc) {
391*60f0b76eSHong Zhang     ierr = MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
392*60f0b76eSHong Zhang     ierr = MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
393*60f0b76eSHong Zhang     ierr = MatSetOption(BB,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
394*60f0b76eSHong Zhang   }
395*60f0b76eSHong Zhang   PetscFunctionReturn(0);
396*60f0b76eSHong Zhang }
397