xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/reaction_diffusion.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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  */
279371c9d4SSatish Balay PetscErrorCode RHSFunction(TS ts, PetscReal ftime, Vec U, Vec F, void *ptr) {
2860f0b76eSHong Zhang   AppCtx     *appctx = (AppCtx *)ptr;
2960f0b76eSHong Zhang   DM          da;
3060f0b76eSHong Zhang   PetscInt    i, j, Mx, My, xs, ys, xm, ym;
3160f0b76eSHong Zhang   PetscReal   hx, hy, sx, sy;
3260f0b76eSHong Zhang   PetscScalar uc, uxx, uyy, vc, vxx, vyy;
3360f0b76eSHong Zhang   Field     **u, **f;
3460f0b76eSHong Zhang   Vec         localU;
3560f0b76eSHong Zhang 
3660f0b76eSHong Zhang   PetscFunctionBegin;
379566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
389566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localU));
399566063dSJacob 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));
409371c9d4SSatish Balay   hx = 2.50 / (PetscReal)Mx;
419371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
429371c9d4SSatish Balay   hy = 2.50 / (PetscReal)My;
439371c9d4SSatish Balay   sy = 1.0 / (hy * hy);
4460f0b76eSHong Zhang 
4560f0b76eSHong Zhang   /*
4660f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
4760f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
4860f0b76eSHong Zhang      By placing code between these two statements, computations can be
4960f0b76eSHong Zhang      done while messages are in transition.
5060f0b76eSHong Zhang   */
519566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
529566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
5360f0b76eSHong Zhang 
5460f0b76eSHong Zhang   /*
5560f0b76eSHong Zhang      Get pointers to vector data
5660f0b76eSHong Zhang   */
579566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localU, &u));
589566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
5960f0b76eSHong Zhang 
6060f0b76eSHong Zhang   /*
6160f0b76eSHong Zhang      Get local grid boundaries
6260f0b76eSHong Zhang   */
639566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
6460f0b76eSHong Zhang 
6560f0b76eSHong Zhang   /*
6660f0b76eSHong Zhang      Compute function over the locally owned part of the grid
6760f0b76eSHong Zhang   */
6860f0b76eSHong Zhang   for (j = ys; j < ys + ym; j++) {
6960f0b76eSHong Zhang     for (i = xs; i < xs + xm; i++) {
7060f0b76eSHong Zhang       uc        = u[j][i].u;
7160f0b76eSHong Zhang       uxx       = (-2.0 * uc + u[j][i - 1].u + u[j][i + 1].u) * sx;
7260f0b76eSHong Zhang       uyy       = (-2.0 * uc + u[j - 1][i].u + u[j + 1][i].u) * sy;
7360f0b76eSHong Zhang       vc        = u[j][i].v;
7460f0b76eSHong Zhang       vxx       = (-2.0 * vc + u[j][i - 1].v + u[j][i + 1].v) * sx;
7560f0b76eSHong Zhang       vyy       = (-2.0 * vc + u[j - 1][i].v + u[j + 1][i].v) * sy;
7660f0b76eSHong Zhang       f[j][i].u = appctx->D1 * (uxx + uyy) - uc * vc * vc + appctx->gamma * (1.0 - uc);
7760f0b76eSHong Zhang       f[j][i].v = appctx->D2 * (vxx + vyy) + uc * vc * vc - (appctx->gamma + appctx->kappa) * vc;
7860f0b76eSHong Zhang     }
7960f0b76eSHong Zhang   }
809566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(16.0 * xm * ym));
8160f0b76eSHong Zhang 
8260f0b76eSHong Zhang   /*
8360f0b76eSHong Zhang      Restore vectors
8460f0b76eSHong Zhang   */
859566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
869566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
879566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localU));
8860f0b76eSHong Zhang   PetscFunctionReturn(0);
8960f0b76eSHong Zhang }
9060f0b76eSHong Zhang 
919371c9d4SSatish Balay PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat BB, void *ctx) {
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;
1039566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1049566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localU));
1059566063dSJacob 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 
1079371c9d4SSatish Balay   hx = 2.50 / (PetscReal)Mx;
1089371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
1099371c9d4SSatish Balay   hy = 2.50 / (PetscReal)My;
1109371c9d4SSatish Balay   sy = 1.0 / (hy * hy);
11160f0b76eSHong Zhang 
11260f0b76eSHong Zhang   /*
11360f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
11460f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
11560f0b76eSHong Zhang      By placing code between these two statements, computations can be
11660f0b76eSHong Zhang      done while messages are in transition.
11760f0b76eSHong Zhang   */
1189566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
1199566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
12060f0b76eSHong Zhang 
12160f0b76eSHong Zhang   /*
12260f0b76eSHong Zhang      Get pointers to vector data
12360f0b76eSHong Zhang   */
1249566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localU, &u));
12560f0b76eSHong Zhang 
12660f0b76eSHong Zhang   /*
12760f0b76eSHong Zhang      Get local grid boundaries
12860f0b76eSHong Zhang   */
1299566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
13060f0b76eSHong Zhang 
13160f0b76eSHong Zhang   stencil[0].k = 0;
13260f0b76eSHong Zhang   stencil[1].k = 0;
13360f0b76eSHong Zhang   stencil[2].k = 0;
13460f0b76eSHong Zhang   stencil[3].k = 0;
13560f0b76eSHong Zhang   stencil[4].k = 0;
13660f0b76eSHong Zhang   stencil[5].k = 0;
13760f0b76eSHong Zhang   rowstencil.k = 0;
13860f0b76eSHong Zhang   /*
13960f0b76eSHong Zhang      Compute function over the locally owned part of the grid
14060f0b76eSHong Zhang   */
14160f0b76eSHong Zhang   for (j = ys; j < ys + ym; j++) {
14260f0b76eSHong Zhang     stencil[0].j = j - 1;
14360f0b76eSHong Zhang     stencil[1].j = j + 1;
14460f0b76eSHong Zhang     stencil[2].j = j;
14560f0b76eSHong Zhang     stencil[3].j = j;
14660f0b76eSHong Zhang     stencil[4].j = j;
14760f0b76eSHong Zhang     stencil[5].j = j;
1489371c9d4SSatish Balay     rowstencil.k = 0;
1499371c9d4SSatish Balay     rowstencil.j = j;
15060f0b76eSHong Zhang     for (i = xs; i < xs + xm; i++) {
15160f0b76eSHong Zhang       uc = u[j][i].u;
15260f0b76eSHong Zhang       vc = u[j][i].v;
15360f0b76eSHong Zhang 
15460f0b76eSHong Zhang       /*      uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
15560f0b76eSHong Zhang       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
15660f0b76eSHong Zhang 
15760f0b76eSHong Zhang       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
15860f0b76eSHong Zhang       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
15960f0b76eSHong Zhang        f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
16060f0b76eSHong Zhang 
1619371c9d4SSatish Balay       stencil[0].i = i;
1629371c9d4SSatish Balay       stencil[0].c = 0;
1639371c9d4SSatish Balay       entries[0]   = appctx->D1 * sy;
1649371c9d4SSatish Balay       stencil[1].i = i;
1659371c9d4SSatish Balay       stencil[1].c = 0;
1669371c9d4SSatish Balay       entries[1]   = appctx->D1 * sy;
1679371c9d4SSatish Balay       stencil[2].i = i - 1;
1689371c9d4SSatish Balay       stencil[2].c = 0;
1699371c9d4SSatish Balay       entries[2]   = appctx->D1 * sx;
1709371c9d4SSatish Balay       stencil[3].i = i + 1;
1719371c9d4SSatish Balay       stencil[3].c = 0;
1729371c9d4SSatish Balay       entries[3]   = appctx->D1 * sx;
1739371c9d4SSatish Balay       stencil[4].i = i;
1749371c9d4SSatish Balay       stencil[4].c = 0;
1759371c9d4SSatish Balay       entries[4]   = -2.0 * appctx->D1 * (sx + sy) - vc * vc - appctx->gamma;
1769371c9d4SSatish Balay       stencil[5].i = i;
1779371c9d4SSatish Balay       stencil[5].c = 1;
1789371c9d4SSatish Balay       entries[5]   = -2.0 * uc * vc;
1799371c9d4SSatish Balay       rowstencil.i = i;
1809371c9d4SSatish Balay       rowstencil.c = 0;
18160f0b76eSHong Zhang 
1829566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
183*48a46eb9SPierre Jolivet       if (appctx->aijpc) PetscCall(MatSetValuesStencil(BB, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
1849371c9d4SSatish Balay       stencil[0].c = 1;
1859371c9d4SSatish Balay       entries[0]   = appctx->D2 * sy;
1869371c9d4SSatish Balay       stencil[1].c = 1;
1879371c9d4SSatish Balay       entries[1]   = appctx->D2 * sy;
1889371c9d4SSatish Balay       stencil[2].c = 1;
1899371c9d4SSatish Balay       entries[2]   = appctx->D2 * sx;
1909371c9d4SSatish Balay       stencil[3].c = 1;
1919371c9d4SSatish Balay       entries[3]   = appctx->D2 * sx;
1929371c9d4SSatish Balay       stencil[4].c = 1;
1939371c9d4SSatish Balay       entries[4]   = -2.0 * appctx->D2 * (sx + sy) + 2.0 * uc * vc - appctx->gamma - appctx->kappa;
1949371c9d4SSatish Balay       stencil[5].c = 0;
1959371c9d4SSatish Balay       entries[5]   = vc * vc;
19660f0b76eSHong Zhang       rowstencil.c = 1;
19760f0b76eSHong Zhang 
1989566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
199*48a46eb9SPierre Jolivet       if (appctx->aijpc) PetscCall(MatSetValuesStencil(BB, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
20060f0b76eSHong Zhang       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
20160f0b76eSHong Zhang     }
20260f0b76eSHong Zhang   }
20360f0b76eSHong Zhang 
20460f0b76eSHong Zhang   /*
20560f0b76eSHong Zhang      Restore vectors
20660f0b76eSHong Zhang   */
2079566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(19.0 * xm * ym));
2089566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
2099566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localU));
2109566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2119566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2129566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
21360f0b76eSHong Zhang   if (appctx->aijpc) {
2149566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(BB, MAT_FINAL_ASSEMBLY));
2159566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(BB, MAT_FINAL_ASSEMBLY));
2169566063dSJacob Faibussowitsch     PetscCall(MatSetOption(BB, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
21760f0b76eSHong Zhang   }
21860f0b76eSHong Zhang   PetscFunctionReturn(0);
21960f0b76eSHong Zhang }
22060f0b76eSHong Zhang 
22160f0b76eSHong Zhang /*
22260f0b76eSHong Zhang    IFunction - Evaluates implicit nonlinear function, xdot - F(x).
22360f0b76eSHong Zhang 
22460f0b76eSHong Zhang    Input Parameters:
22560f0b76eSHong Zhang .  ts - the TS context
22660f0b76eSHong Zhang .  U - input vector
22760f0b76eSHong Zhang .  Udot - input vector
22860f0b76eSHong Zhang .  ptr - optional user-defined context, as set by TSSetRHSFunction()
22960f0b76eSHong Zhang 
23060f0b76eSHong Zhang    Output Parameter:
23160f0b76eSHong Zhang .  F - function vector
23260f0b76eSHong Zhang  */
2339371c9d4SSatish Balay PetscErrorCode IFunction(TS ts, PetscReal ftime, Vec U, Vec Udot, Vec F, void *ptr) {
23460f0b76eSHong Zhang   AppCtx     *appctx = (AppCtx *)ptr;
23560f0b76eSHong Zhang   DM          da;
23660f0b76eSHong Zhang   PetscInt    i, j, Mx, My, xs, ys, xm, ym;
23760f0b76eSHong Zhang   PetscReal   hx, hy, sx, sy;
23860f0b76eSHong Zhang   PetscScalar uc, uxx, uyy, vc, vxx, vyy;
23960f0b76eSHong Zhang   Field     **u, **f, **udot;
24060f0b76eSHong Zhang   Vec         localU;
24160f0b76eSHong Zhang 
24260f0b76eSHong Zhang   PetscFunctionBegin;
2439566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2449566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localU));
2459566063dSJacob 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));
2469371c9d4SSatish Balay   hx = 2.50 / (PetscReal)Mx;
2479371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
2489371c9d4SSatish Balay   hy = 2.50 / (PetscReal)My;
2499371c9d4SSatish Balay   sy = 1.0 / (hy * hy);
25060f0b76eSHong Zhang 
25160f0b76eSHong Zhang   /*
25260f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
25360f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
25460f0b76eSHong Zhang      By placing code between these two statements, computations can be
25560f0b76eSHong Zhang      done while messages are in transition.
25660f0b76eSHong Zhang   */
2579566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
2589566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
25960f0b76eSHong Zhang 
26060f0b76eSHong Zhang   /*
26160f0b76eSHong Zhang      Get pointers to vector data
26260f0b76eSHong Zhang   */
2639566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localU, &u));
2649566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
2659566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, Udot, &udot));
26660f0b76eSHong Zhang 
26760f0b76eSHong Zhang   /*
26860f0b76eSHong Zhang      Get local grid boundaries
26960f0b76eSHong Zhang   */
2709566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
27160f0b76eSHong Zhang 
27260f0b76eSHong Zhang   /*
27360f0b76eSHong Zhang      Compute function over the locally owned part of the grid
27460f0b76eSHong Zhang   */
27560f0b76eSHong Zhang   for (j = ys; j < ys + ym; j++) {
27660f0b76eSHong Zhang     for (i = xs; i < xs + xm; i++) {
27760f0b76eSHong Zhang       uc        = u[j][i].u;
27860f0b76eSHong Zhang       uxx       = (-2.0 * uc + u[j][i - 1].u + u[j][i + 1].u) * sx;
27960f0b76eSHong Zhang       uyy       = (-2.0 * uc + u[j - 1][i].u + u[j + 1][i].u) * sy;
28060f0b76eSHong Zhang       vc        = u[j][i].v;
28160f0b76eSHong Zhang       vxx       = (-2.0 * vc + u[j][i - 1].v + u[j][i + 1].v) * sx;
28260f0b76eSHong Zhang       vyy       = (-2.0 * vc + u[j - 1][i].v + u[j + 1][i].v) * sy;
28360f0b76eSHong Zhang       f[j][i].u = udot[j][i].u - (appctx->D1 * (uxx + uyy) - uc * vc * vc + appctx->gamma * (1.0 - uc));
28460f0b76eSHong Zhang       f[j][i].v = udot[j][i].v - (appctx->D2 * (vxx + vyy) + uc * vc * vc - (appctx->gamma + appctx->kappa) * vc);
28560f0b76eSHong Zhang     }
28660f0b76eSHong Zhang   }
2879566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(16.0 * xm * ym));
28860f0b76eSHong Zhang 
28960f0b76eSHong Zhang   /*
29060f0b76eSHong Zhang      Restore vectors
29160f0b76eSHong Zhang   */
2929566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
2939566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
2949566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, Udot, &udot));
2959566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localU));
29660f0b76eSHong Zhang   PetscFunctionReturn(0);
29760f0b76eSHong Zhang }
29860f0b76eSHong Zhang 
2999371c9d4SSatish Balay PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat BB, void *ctx) {
30060f0b76eSHong Zhang   AppCtx     *appctx = (AppCtx *)ctx; /* user-defined application context */
30160f0b76eSHong Zhang   DM          da;
30260f0b76eSHong Zhang   PetscInt    i, j, Mx, My, xs, ys, xm, ym;
30360f0b76eSHong Zhang   PetscReal   hx, hy, sx, sy;
30460f0b76eSHong Zhang   PetscScalar uc, vc;
30560f0b76eSHong Zhang   Field     **u;
30660f0b76eSHong Zhang   Vec         localU;
30760f0b76eSHong Zhang   MatStencil  stencil[6], rowstencil;
30860f0b76eSHong Zhang   PetscScalar entries[6];
30960f0b76eSHong Zhang 
31060f0b76eSHong Zhang   PetscFunctionBegin;
3119566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
3129566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localU));
3139566063dSJacob 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));
31460f0b76eSHong Zhang 
3159371c9d4SSatish Balay   hx = 2.50 / (PetscReal)Mx;
3169371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
3179371c9d4SSatish Balay   hy = 2.50 / (PetscReal)My;
3189371c9d4SSatish Balay   sy = 1.0 / (hy * hy);
31960f0b76eSHong Zhang 
32060f0b76eSHong Zhang   /*
32160f0b76eSHong Zhang      Scatter ghost points to local vector,using the 2-step process
32260f0b76eSHong Zhang         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
32360f0b76eSHong Zhang      By placing code between these two statements, computations can be
32460f0b76eSHong Zhang      done while messages are in transition.
32560f0b76eSHong Zhang   */
3269566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
3279566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
32860f0b76eSHong Zhang 
32960f0b76eSHong Zhang   /*
33060f0b76eSHong Zhang      Get pointers to vector data
33160f0b76eSHong Zhang   */
3329566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localU, &u));
33360f0b76eSHong Zhang 
33460f0b76eSHong Zhang   /*
33560f0b76eSHong Zhang      Get local grid boundaries
33660f0b76eSHong Zhang   */
3379566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
33860f0b76eSHong Zhang 
33960f0b76eSHong Zhang   stencil[0].k = 0;
34060f0b76eSHong Zhang   stencil[1].k = 0;
34160f0b76eSHong Zhang   stencil[2].k = 0;
34260f0b76eSHong Zhang   stencil[3].k = 0;
34360f0b76eSHong Zhang   stencil[4].k = 0;
34460f0b76eSHong Zhang   stencil[5].k = 0;
34560f0b76eSHong Zhang   rowstencil.k = 0;
34660f0b76eSHong Zhang   /*
34760f0b76eSHong Zhang      Compute function over the locally owned part of the grid
34860f0b76eSHong Zhang   */
34960f0b76eSHong Zhang   for (j = ys; j < ys + ym; j++) {
35060f0b76eSHong Zhang     stencil[0].j = j - 1;
35160f0b76eSHong Zhang     stencil[1].j = j + 1;
35260f0b76eSHong Zhang     stencil[2].j = j;
35360f0b76eSHong Zhang     stencil[3].j = j;
35460f0b76eSHong Zhang     stencil[4].j = j;
35560f0b76eSHong Zhang     stencil[5].j = j;
3569371c9d4SSatish Balay     rowstencil.k = 0;
3579371c9d4SSatish Balay     rowstencil.j = j;
35860f0b76eSHong Zhang     for (i = xs; i < xs + xm; i++) {
35960f0b76eSHong Zhang       uc = u[j][i].u;
36060f0b76eSHong Zhang       vc = u[j][i].v;
36160f0b76eSHong Zhang 
36260f0b76eSHong Zhang       /*      uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
36360f0b76eSHong Zhang       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
36460f0b76eSHong Zhang 
36560f0b76eSHong Zhang       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
36660f0b76eSHong Zhang       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
36760f0b76eSHong Zhang        f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
36860f0b76eSHong Zhang 
3699371c9d4SSatish Balay       stencil[0].i = i;
3709371c9d4SSatish Balay       stencil[0].c = 0;
3719371c9d4SSatish Balay       entries[0]   = -appctx->D1 * sy;
3729371c9d4SSatish Balay       stencil[1].i = i;
3739371c9d4SSatish Balay       stencil[1].c = 0;
3749371c9d4SSatish Balay       entries[1]   = -appctx->D1 * sy;
3759371c9d4SSatish Balay       stencil[2].i = i - 1;
3769371c9d4SSatish Balay       stencil[2].c = 0;
3779371c9d4SSatish Balay       entries[2]   = -appctx->D1 * sx;
3789371c9d4SSatish Balay       stencil[3].i = i + 1;
3799371c9d4SSatish Balay       stencil[3].c = 0;
3809371c9d4SSatish Balay       entries[3]   = -appctx->D1 * sx;
3819371c9d4SSatish Balay       stencil[4].i = i;
3829371c9d4SSatish Balay       stencil[4].c = 0;
3839371c9d4SSatish Balay       entries[4]   = 2.0 * appctx->D1 * (sx + sy) + vc * vc + appctx->gamma + a;
3849371c9d4SSatish Balay       stencil[5].i = i;
3859371c9d4SSatish Balay       stencil[5].c = 1;
3869371c9d4SSatish Balay       entries[5]   = 2.0 * uc * vc;
3879371c9d4SSatish Balay       rowstencil.i = i;
3889371c9d4SSatish Balay       rowstencil.c = 0;
38960f0b76eSHong Zhang 
3909566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
391*48a46eb9SPierre Jolivet       if (appctx->aijpc) PetscCall(MatSetValuesStencil(BB, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
3929371c9d4SSatish Balay       stencil[0].c = 1;
3939371c9d4SSatish Balay       entries[0]   = -appctx->D2 * sy;
3949371c9d4SSatish Balay       stencil[1].c = 1;
3959371c9d4SSatish Balay       entries[1]   = -appctx->D2 * sy;
3969371c9d4SSatish Balay       stencil[2].c = 1;
3979371c9d4SSatish Balay       entries[2]   = -appctx->D2 * sx;
3989371c9d4SSatish Balay       stencil[3].c = 1;
3999371c9d4SSatish Balay       entries[3]   = -appctx->D2 * sx;
4009371c9d4SSatish Balay       stencil[4].c = 1;
4019371c9d4SSatish Balay       entries[4]   = 2.0 * appctx->D2 * (sx + sy) - 2.0 * uc * vc + appctx->gamma + appctx->kappa + a;
4029371c9d4SSatish Balay       stencil[5].c = 0;
4039371c9d4SSatish Balay       entries[5]   = -vc * vc;
40460f0b76eSHong Zhang       rowstencil.c = 1;
40560f0b76eSHong Zhang 
4069566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
407*48a46eb9SPierre Jolivet       if (appctx->aijpc) PetscCall(MatSetValuesStencil(BB, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
40860f0b76eSHong Zhang       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
40960f0b76eSHong Zhang     }
41060f0b76eSHong Zhang   }
41160f0b76eSHong Zhang 
41260f0b76eSHong Zhang   /*
41360f0b76eSHong Zhang      Restore vectors
41460f0b76eSHong Zhang   */
4159566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(19.0 * xm * ym));
4169566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
4179566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localU));
4189566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
4209566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
42160f0b76eSHong Zhang   if (appctx->aijpc) {
4229566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(BB, MAT_FINAL_ASSEMBLY));
4239566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(BB, MAT_FINAL_ASSEMBLY));
4249566063dSJacob Faibussowitsch     PetscCall(MatSetOption(BB, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
42560f0b76eSHong Zhang   }
42660f0b76eSHong Zhang   PetscFunctionReturn(0);
42760f0b76eSHong Zhang }
428