xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/reaction_diffusion.c (revision 9566063d113dddea24716c546802770db7481bc0)
160f0b76eSHong Zhang #include "reaction_diffusion.h"
260f0b76eSHong Zhang #include <petscdm.h>
360f0b76eSHong Zhang #include <petscdmda.h>
460f0b76eSHong Zhang 
560f0b76eSHong Zhang /*F
660f0b76eSHong Zhang      This example is taken from the book, Numerical Solution of Time-Dependent Advection-Diffusion-Reaction Equations by
760f0b76eSHong Zhang       W. Hundsdorf and J.G. Verwer,  Page 21, Pattern Formation with Reaction-Diffusion Equations
860f0b76eSHong Zhang \begin{eqnarray*}
960f0b76eSHong Zhang         u_t = D_1 (u_{xx} + u_{yy})  - u*v^2 + \gamma(1 -u)           \\
1060f0b76eSHong Zhang         v_t = D_2 (v_{xx} + v_{yy})  + u*v^2 - (\gamma + \kappa)v
1160f0b76eSHong Zhang \end{eqnarray*}
1260f0b76eSHong Zhang     Unlike in the book this uses periodic boundary conditions instead of Neumann
1360f0b76eSHong Zhang     (since they are easier for finite differences).
1460f0b76eSHong Zhang F*/
1560f0b76eSHong Zhang 
1660f0b76eSHong Zhang /*
1760f0b76eSHong Zhang    RHSFunction - Evaluates nonlinear function, F(x).
1860f0b76eSHong Zhang 
1960f0b76eSHong Zhang    Input Parameters:
2060f0b76eSHong Zhang .  ts - the TS context
2160f0b76eSHong Zhang .  X - input vector
2260f0b76eSHong Zhang .  ptr - optional user-defined context, as set by TSSetRHSFunction()
2360f0b76eSHong Zhang 
2460f0b76eSHong Zhang    Output Parameter:
2560f0b76eSHong Zhang .  F - function vector
2660f0b76eSHong Zhang  */
2760f0b76eSHong Zhang PetscErrorCode RHSFunction(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr)
2860f0b76eSHong Zhang {
2960f0b76eSHong Zhang   AppCtx         *appctx = (AppCtx*)ptr;
3060f0b76eSHong Zhang   DM             da;
3160f0b76eSHong Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
3260f0b76eSHong Zhang   PetscReal      hx,hy,sx,sy;
3360f0b76eSHong Zhang   PetscScalar    uc,uxx,uyy,vc,vxx,vyy;
3460f0b76eSHong Zhang   Field          **u,**f;
3560f0b76eSHong Zhang   Vec            localU;
3660f0b76eSHong Zhang 
3760f0b76eSHong Zhang   PetscFunctionBegin;
38*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
39*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
40*9566063dSJacob Faibussowitsch   PetscCall(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));
4160f0b76eSHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
4260f0b76eSHong Zhang   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
4360f0b76eSHong Zhang 
4460f0b76eSHong Zhang   /*
4560f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
4660f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
4760f0b76eSHong Zhang      By placing code between these two statements, computations can be
4860f0b76eSHong Zhang      done while messages are in transition.
4960f0b76eSHong Zhang   */
50*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
51*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
5260f0b76eSHong Zhang 
5360f0b76eSHong Zhang   /*
5460f0b76eSHong Zhang      Get pointers to vector data
5560f0b76eSHong Zhang   */
56*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
57*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,F,&f));
5860f0b76eSHong Zhang 
5960f0b76eSHong Zhang   /*
6060f0b76eSHong Zhang      Get local grid boundaries
6160f0b76eSHong Zhang   */
62*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
6360f0b76eSHong Zhang 
6460f0b76eSHong Zhang   /*
6560f0b76eSHong Zhang      Compute function over the locally owned part of the grid
6660f0b76eSHong Zhang   */
6760f0b76eSHong Zhang   for (j=ys; j<ys+ym; j++) {
6860f0b76eSHong Zhang     for (i=xs; i<xs+xm; i++) {
6960f0b76eSHong Zhang       uc        = u[j][i].u;
7060f0b76eSHong Zhang       uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
7160f0b76eSHong Zhang       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
7260f0b76eSHong Zhang       vc        = u[j][i].v;
7360f0b76eSHong Zhang       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
7460f0b76eSHong Zhang       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
7560f0b76eSHong Zhang       f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);
7660f0b76eSHong Zhang       f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc;
7760f0b76eSHong Zhang     }
7860f0b76eSHong Zhang   }
79*9566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(16.0*xm*ym));
8060f0b76eSHong Zhang 
8160f0b76eSHong Zhang   /*
8260f0b76eSHong Zhang      Restore vectors
8360f0b76eSHong Zhang   */
84*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
85*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,F,&f));
86*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
8760f0b76eSHong Zhang   PetscFunctionReturn(0);
8860f0b76eSHong Zhang }
8960f0b76eSHong Zhang 
9060f0b76eSHong Zhang PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat BB,void *ctx)
9160f0b76eSHong Zhang {
9260f0b76eSHong Zhang   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
9360f0b76eSHong Zhang   DM             da;
9460f0b76eSHong Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
9560f0b76eSHong Zhang   PetscReal      hx,hy,sx,sy;
9660f0b76eSHong Zhang   PetscScalar    uc,vc;
9760f0b76eSHong Zhang   Field          **u;
9860f0b76eSHong Zhang   Vec            localU;
9960f0b76eSHong Zhang   MatStencil     stencil[6],rowstencil;
10060f0b76eSHong Zhang   PetscScalar    entries[6];
10160f0b76eSHong Zhang 
10260f0b76eSHong Zhang   PetscFunctionBegin;
103*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
104*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
105*9566063dSJacob Faibussowitsch   PetscCall(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));
10660f0b76eSHong Zhang 
10760f0b76eSHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
10860f0b76eSHong Zhang   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
10960f0b76eSHong Zhang 
11060f0b76eSHong Zhang   /*
11160f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
11260f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
11360f0b76eSHong Zhang      By placing code between these two statements, computations can be
11460f0b76eSHong Zhang      done while messages are in transition.
11560f0b76eSHong Zhang   */
116*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
117*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
11860f0b76eSHong Zhang 
11960f0b76eSHong Zhang   /*
12060f0b76eSHong Zhang      Get pointers to vector data
12160f0b76eSHong Zhang   */
122*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
12360f0b76eSHong Zhang 
12460f0b76eSHong Zhang   /*
12560f0b76eSHong Zhang      Get local grid boundaries
12660f0b76eSHong Zhang   */
127*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
12860f0b76eSHong Zhang 
12960f0b76eSHong Zhang   stencil[0].k = 0;
13060f0b76eSHong Zhang   stencil[1].k = 0;
13160f0b76eSHong Zhang   stencil[2].k = 0;
13260f0b76eSHong Zhang   stencil[3].k = 0;
13360f0b76eSHong Zhang   stencil[4].k = 0;
13460f0b76eSHong Zhang   stencil[5].k = 0;
13560f0b76eSHong Zhang   rowstencil.k = 0;
13660f0b76eSHong Zhang   /*
13760f0b76eSHong Zhang      Compute function over the locally owned part of the grid
13860f0b76eSHong Zhang   */
13960f0b76eSHong Zhang   for (j=ys; j<ys+ym; j++) {
14060f0b76eSHong Zhang     stencil[0].j = j-1;
14160f0b76eSHong Zhang     stencil[1].j = j+1;
14260f0b76eSHong Zhang     stencil[2].j = j;
14360f0b76eSHong Zhang     stencil[3].j = j;
14460f0b76eSHong Zhang     stencil[4].j = j;
14560f0b76eSHong Zhang     stencil[5].j = j;
14660f0b76eSHong Zhang     rowstencil.k = 0; rowstencil.j = j;
14760f0b76eSHong Zhang     for (i=xs; i<xs+xm; i++) {
14860f0b76eSHong Zhang       uc = u[j][i].u;
14960f0b76eSHong Zhang       vc = u[j][i].v;
15060f0b76eSHong Zhang 
15160f0b76eSHong Zhang       /*      uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
15260f0b76eSHong Zhang       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
15360f0b76eSHong Zhang 
15460f0b76eSHong Zhang       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
15560f0b76eSHong Zhang       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
15660f0b76eSHong Zhang        f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
15760f0b76eSHong Zhang 
15860f0b76eSHong Zhang       stencil[0].i = i; stencil[0].c = 0; entries[0] = appctx->D1*sy;
15960f0b76eSHong Zhang       stencil[1].i = i; stencil[1].c = 0; entries[1] = appctx->D1*sy;
16060f0b76eSHong Zhang       stencil[2].i = i-1; stencil[2].c = 0; entries[2] = appctx->D1*sx;
16160f0b76eSHong Zhang       stencil[3].i = i+1; stencil[3].c = 0; entries[3] = appctx->D1*sx;
16260f0b76eSHong Zhang       stencil[4].i = i; stencil[4].c = 0; entries[4] = -2.0*appctx->D1*(sx + sy) - vc*vc - appctx->gamma;
16360f0b76eSHong Zhang       stencil[5].i = i; stencil[5].c = 1; entries[5] = -2.0*uc*vc;
16460f0b76eSHong Zhang       rowstencil.i = i; rowstencil.c = 0;
16560f0b76eSHong Zhang 
166*9566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
16760f0b76eSHong Zhang       if (appctx->aijpc) {
168*9566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
16960f0b76eSHong Zhang       }
17060f0b76eSHong Zhang       stencil[0].c = 1; entries[0] = appctx->D2*sy;
17160f0b76eSHong Zhang       stencil[1].c = 1; entries[1] = appctx->D2*sy;
17260f0b76eSHong Zhang       stencil[2].c = 1; entries[2] = appctx->D2*sx;
17360f0b76eSHong Zhang       stencil[3].c = 1; entries[3] = appctx->D2*sx;
17460f0b76eSHong Zhang       stencil[4].c = 1; entries[4] = -2.0*appctx->D2*(sx + sy) + 2.0*uc*vc - appctx->gamma - appctx->kappa;
17560f0b76eSHong Zhang       stencil[5].c = 0; entries[5] = vc*vc;
17660f0b76eSHong Zhang       rowstencil.c = 1;
17760f0b76eSHong Zhang 
178*9566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
17960f0b76eSHong Zhang       if (appctx->aijpc) {
180*9566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
18160f0b76eSHong Zhang       }
18260f0b76eSHong Zhang       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
18360f0b76eSHong Zhang     }
18460f0b76eSHong Zhang   }
18560f0b76eSHong Zhang 
18660f0b76eSHong Zhang   /*
18760f0b76eSHong Zhang      Restore vectors
18860f0b76eSHong Zhang   */
189*9566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(19.0*xm*ym));
190*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
191*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
192*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
193*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
194*9566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
19560f0b76eSHong Zhang   if (appctx->aijpc) {
196*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY));
197*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY));
198*9566063dSJacob Faibussowitsch     PetscCall(MatSetOption(BB,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
19960f0b76eSHong Zhang   }
20060f0b76eSHong Zhang   PetscFunctionReturn(0);
20160f0b76eSHong Zhang }
20260f0b76eSHong Zhang 
20360f0b76eSHong Zhang /*
20460f0b76eSHong Zhang    IFunction - Evaluates implicit nonlinear function, xdot - F(x).
20560f0b76eSHong Zhang 
20660f0b76eSHong Zhang    Input Parameters:
20760f0b76eSHong Zhang .  ts - the TS context
20860f0b76eSHong Zhang .  U - input vector
20960f0b76eSHong Zhang .  Udot - input vector
21060f0b76eSHong Zhang .  ptr - optional user-defined context, as set by TSSetRHSFunction()
21160f0b76eSHong Zhang 
21260f0b76eSHong Zhang    Output Parameter:
21360f0b76eSHong Zhang .  F - function vector
21460f0b76eSHong Zhang  */
21560f0b76eSHong Zhang PetscErrorCode IFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
21660f0b76eSHong Zhang {
21760f0b76eSHong Zhang   AppCtx         *appctx = (AppCtx*)ptr;
21860f0b76eSHong Zhang   DM             da;
21960f0b76eSHong Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
22060f0b76eSHong Zhang   PetscReal      hx,hy,sx,sy;
22160f0b76eSHong Zhang   PetscScalar    uc,uxx,uyy,vc,vxx,vyy;
22260f0b76eSHong Zhang   Field          **u,**f,**udot;
22360f0b76eSHong Zhang   Vec            localU;
22460f0b76eSHong Zhang 
22560f0b76eSHong Zhang   PetscFunctionBegin;
226*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
227*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
228*9566063dSJacob Faibussowitsch   PetscCall(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));
22960f0b76eSHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
23060f0b76eSHong Zhang   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
23160f0b76eSHong Zhang 
23260f0b76eSHong Zhang   /*
23360f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
23460f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
23560f0b76eSHong Zhang      By placing code between these two statements, computations can be
23660f0b76eSHong Zhang      done while messages are in transition.
23760f0b76eSHong Zhang   */
238*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
239*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
24060f0b76eSHong Zhang 
24160f0b76eSHong Zhang   /*
24260f0b76eSHong Zhang      Get pointers to vector data
24360f0b76eSHong Zhang   */
244*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
245*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,F,&f));
246*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,Udot,&udot));
24760f0b76eSHong Zhang 
24860f0b76eSHong Zhang   /*
24960f0b76eSHong Zhang      Get local grid boundaries
25060f0b76eSHong Zhang   */
251*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
25260f0b76eSHong Zhang 
25360f0b76eSHong Zhang   /*
25460f0b76eSHong Zhang      Compute function over the locally owned part of the grid
25560f0b76eSHong Zhang   */
25660f0b76eSHong Zhang   for (j=ys; j<ys+ym; j++) {
25760f0b76eSHong Zhang     for (i=xs; i<xs+xm; i++) {
25860f0b76eSHong Zhang       uc        = u[j][i].u;
25960f0b76eSHong Zhang       uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
26060f0b76eSHong Zhang       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
26160f0b76eSHong Zhang       vc        = u[j][i].v;
26260f0b76eSHong Zhang       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
26360f0b76eSHong Zhang       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
26460f0b76eSHong Zhang       f[j][i].u = udot[j][i].u - ( appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc));
26560f0b76eSHong Zhang       f[j][i].v = udot[j][i].v - ( appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc);
26660f0b76eSHong Zhang     }
26760f0b76eSHong Zhang   }
268*9566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(16.0*xm*ym));
26960f0b76eSHong Zhang 
27060f0b76eSHong Zhang   /*
27160f0b76eSHong Zhang      Restore vectors
27260f0b76eSHong Zhang   */
273*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
274*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,F,&f));
275*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,Udot,&udot));
276*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
27760f0b76eSHong Zhang   PetscFunctionReturn(0);
27860f0b76eSHong Zhang }
27960f0b76eSHong Zhang 
28060f0b76eSHong Zhang PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat BB,void *ctx)
28160f0b76eSHong Zhang {
28260f0b76eSHong Zhang   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
28360f0b76eSHong Zhang   DM             da;
28460f0b76eSHong Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
28560f0b76eSHong Zhang   PetscReal      hx,hy,sx,sy;
28660f0b76eSHong Zhang   PetscScalar    uc,vc;
28760f0b76eSHong Zhang   Field          **u;
28860f0b76eSHong Zhang   Vec            localU;
28960f0b76eSHong Zhang   MatStencil     stencil[6],rowstencil;
29060f0b76eSHong Zhang   PetscScalar    entries[6];
29160f0b76eSHong Zhang 
29260f0b76eSHong Zhang   PetscFunctionBegin;
293*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
294*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
295*9566063dSJacob Faibussowitsch   PetscCall(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));
29660f0b76eSHong Zhang 
29760f0b76eSHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
29860f0b76eSHong Zhang   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
29960f0b76eSHong Zhang 
30060f0b76eSHong Zhang   /*
30160f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
30260f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
30360f0b76eSHong Zhang      By placing code between these two statements, computations can be
30460f0b76eSHong Zhang      done while messages are in transition.
30560f0b76eSHong Zhang   */
306*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
307*9566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
30860f0b76eSHong Zhang 
30960f0b76eSHong Zhang   /*
31060f0b76eSHong Zhang      Get pointers to vector data
31160f0b76eSHong Zhang   */
312*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
31360f0b76eSHong Zhang 
31460f0b76eSHong Zhang   /*
31560f0b76eSHong Zhang      Get local grid boundaries
31660f0b76eSHong Zhang   */
317*9566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
31860f0b76eSHong Zhang 
31960f0b76eSHong Zhang   stencil[0].k = 0;
32060f0b76eSHong Zhang   stencil[1].k = 0;
32160f0b76eSHong Zhang   stencil[2].k = 0;
32260f0b76eSHong Zhang   stencil[3].k = 0;
32360f0b76eSHong Zhang   stencil[4].k = 0;
32460f0b76eSHong Zhang   stencil[5].k = 0;
32560f0b76eSHong Zhang   rowstencil.k = 0;
32660f0b76eSHong Zhang   /*
32760f0b76eSHong Zhang      Compute function over the locally owned part of the grid
32860f0b76eSHong Zhang   */
32960f0b76eSHong Zhang   for (j=ys; j<ys+ym; j++) {
33060f0b76eSHong Zhang 
33160f0b76eSHong Zhang     stencil[0].j = j-1;
33260f0b76eSHong Zhang     stencil[1].j = j+1;
33360f0b76eSHong Zhang     stencil[2].j = j;
33460f0b76eSHong Zhang     stencil[3].j = j;
33560f0b76eSHong Zhang     stencil[4].j = j;
33660f0b76eSHong Zhang     stencil[5].j = j;
33760f0b76eSHong Zhang     rowstencil.k = 0; rowstencil.j = j;
33860f0b76eSHong Zhang     for (i=xs; i<xs+xm; i++) {
33960f0b76eSHong Zhang       uc = u[j][i].u;
34060f0b76eSHong Zhang       vc = u[j][i].v;
34160f0b76eSHong Zhang 
34260f0b76eSHong Zhang       /*      uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
34360f0b76eSHong Zhang       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
34460f0b76eSHong Zhang 
34560f0b76eSHong Zhang       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
34660f0b76eSHong Zhang       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
34760f0b76eSHong Zhang        f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
34860f0b76eSHong Zhang 
34960f0b76eSHong Zhang       stencil[0].i = i; stencil[0].c = 0; entries[0] = -appctx->D1*sy;
35060f0b76eSHong Zhang       stencil[1].i = i; stencil[1].c = 0; entries[1] = -appctx->D1*sy;
35160f0b76eSHong Zhang       stencil[2].i = i-1; stencil[2].c = 0; entries[2] = -appctx->D1*sx;
35260f0b76eSHong Zhang       stencil[3].i = i+1; stencil[3].c = 0; entries[3] = -appctx->D1*sx;
35360f0b76eSHong Zhang       stencil[4].i = i; stencil[4].c = 0; entries[4] = 2.0*appctx->D1*(sx + sy) + vc*vc + appctx->gamma + a;
35460f0b76eSHong Zhang       stencil[5].i = i; stencil[5].c = 1; entries[5] = 2.0*uc*vc;
35560f0b76eSHong Zhang       rowstencil.i = i; rowstencil.c = 0;
35660f0b76eSHong Zhang 
357*9566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
35860f0b76eSHong Zhang       if (appctx->aijpc) {
359*9566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
36060f0b76eSHong Zhang       }
36160f0b76eSHong Zhang       stencil[0].c = 1; entries[0] = -appctx->D2*sy;
36260f0b76eSHong Zhang       stencil[1].c = 1; entries[1] = -appctx->D2*sy;
36360f0b76eSHong Zhang       stencil[2].c = 1; entries[2] = -appctx->D2*sx;
36460f0b76eSHong Zhang       stencil[3].c = 1; entries[3] = -appctx->D2*sx;
36560f0b76eSHong Zhang       stencil[4].c = 1; entries[4] = 2.0*appctx->D2*(sx + sy) - 2.0*uc*vc + appctx->gamma + appctx->kappa + a;
36660f0b76eSHong Zhang       stencil[5].c = 0; entries[5] = -vc*vc;
36760f0b76eSHong Zhang       rowstencil.c = 1;
36860f0b76eSHong Zhang 
369*9566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
37060f0b76eSHong Zhang       if (appctx->aijpc) {
371*9566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(BB,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
37260f0b76eSHong Zhang       }
37360f0b76eSHong Zhang       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
37460f0b76eSHong Zhang     }
37560f0b76eSHong Zhang   }
37660f0b76eSHong Zhang 
37760f0b76eSHong Zhang   /*
37860f0b76eSHong Zhang      Restore vectors
37960f0b76eSHong Zhang   */
380*9566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(19.0*xm*ym));
381*9566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
382*9566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
383*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
384*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
385*9566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
38660f0b76eSHong Zhang   if (appctx->aijpc) {
387*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY));
388*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY));
389*9566063dSJacob Faibussowitsch     PetscCall(MatSetOption(BB,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
39060f0b76eSHong Zhang   }
39160f0b76eSHong Zhang   PetscFunctionReturn(0);
39260f0b76eSHong Zhang }
393