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 */
RHSFunction(TS ts,PetscReal ftime,Vec U,Vec F,void * ptr)27d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal ftime, Vec U, Vec F, void *ptr)
28d71ae5a4SJacob Faibussowitsch {
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;
389566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
399566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
409566063dSJacob 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));
419371c9d4SSatish Balay hx = 2.50 / (PetscReal)Mx;
429371c9d4SSatish Balay sx = 1.0 / (hx * hx);
439371c9d4SSatish Balay hy = 2.50 / (PetscReal)My;
449371c9d4SSatish Balay sy = 1.0 / (hy * hy);
4560f0b76eSHong Zhang
4660f0b76eSHong Zhang /*
4760f0b76eSHong Zhang Scatter ghost points to local vector,using the 2-step process
4860f0b76eSHong Zhang DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
4960f0b76eSHong Zhang By placing code between these two statements, computations can be
5060f0b76eSHong Zhang done while messages are in transition.
5160f0b76eSHong Zhang */
529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
5460f0b76eSHong Zhang
5560f0b76eSHong Zhang /*
5660f0b76eSHong Zhang Get pointers to vector data
5760f0b76eSHong Zhang */
589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u));
599566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f));
6060f0b76eSHong Zhang
6160f0b76eSHong Zhang /*
6260f0b76eSHong Zhang Get local grid boundaries
6360f0b76eSHong Zhang */
649566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
6560f0b76eSHong Zhang
6660f0b76eSHong Zhang /*
6760f0b76eSHong Zhang Compute function over the locally owned part of the grid
6860f0b76eSHong Zhang */
6960f0b76eSHong Zhang for (j = ys; j < ys + ym; j++) {
7060f0b76eSHong Zhang for (i = xs; i < xs + xm; i++) {
7160f0b76eSHong Zhang uc = u[j][i].u;
7260f0b76eSHong Zhang uxx = (-2.0 * uc + u[j][i - 1].u + u[j][i + 1].u) * sx;
7360f0b76eSHong Zhang uyy = (-2.0 * uc + u[j - 1][i].u + u[j + 1][i].u) * sy;
7460f0b76eSHong Zhang vc = u[j][i].v;
7560f0b76eSHong Zhang vxx = (-2.0 * vc + u[j][i - 1].v + u[j][i + 1].v) * sx;
7660f0b76eSHong Zhang vyy = (-2.0 * vc + u[j - 1][i].v + u[j + 1][i].v) * sy;
7760f0b76eSHong Zhang f[j][i].u = appctx->D1 * (uxx + uyy) - uc * vc * vc + appctx->gamma * (1.0 - uc);
7860f0b76eSHong Zhang f[j][i].v = appctx->D2 * (vxx + vyy) + uc * vc * vc - (appctx->gamma + appctx->kappa) * vc;
7960f0b76eSHong Zhang }
8060f0b76eSHong Zhang }
819566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym));
8260f0b76eSHong Zhang
8360f0b76eSHong Zhang /*
8460f0b76eSHong Zhang Restore vectors
8560f0b76eSHong Zhang */
869566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
879566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f));
889566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9060f0b76eSHong Zhang }
9160f0b76eSHong Zhang
RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat BB,PetscCtx ctx)92*2a8381b2SBarry Smith PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat BB, PetscCtx ctx)
93d71ae5a4SJacob Faibussowitsch {
9460f0b76eSHong Zhang AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
9560f0b76eSHong Zhang DM da;
9660f0b76eSHong Zhang PetscInt i, j, Mx, My, xs, ys, xm, ym;
9760f0b76eSHong Zhang PetscReal hx, hy, sx, sy;
9860f0b76eSHong Zhang PetscScalar uc, vc;
9960f0b76eSHong Zhang Field **u;
10060f0b76eSHong Zhang Vec localU;
10160f0b76eSHong Zhang MatStencil stencil[6], rowstencil;
10260f0b76eSHong Zhang PetscScalar entries[6];
10360f0b76eSHong Zhang
10460f0b76eSHong Zhang PetscFunctionBegin;
1059566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
1069566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
1079566063dSJacob 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));
10860f0b76eSHong Zhang
1099371c9d4SSatish Balay hx = 2.50 / (PetscReal)Mx;
1109371c9d4SSatish Balay sx = 1.0 / (hx * hx);
1119371c9d4SSatish Balay hy = 2.50 / (PetscReal)My;
1129371c9d4SSatish Balay sy = 1.0 / (hy * hy);
11360f0b76eSHong Zhang
11460f0b76eSHong Zhang /*
11560f0b76eSHong Zhang Scatter ghost points to local vector,using the 2-step process
11660f0b76eSHong Zhang DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
11760f0b76eSHong Zhang By placing code between these two statements, computations can be
11860f0b76eSHong Zhang done while messages are in transition.
11960f0b76eSHong Zhang */
1209566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
1219566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
12260f0b76eSHong Zhang
12360f0b76eSHong Zhang /*
12460f0b76eSHong Zhang Get pointers to vector data
12560f0b76eSHong Zhang */
1269566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u));
12760f0b76eSHong Zhang
12860f0b76eSHong Zhang /*
12960f0b76eSHong Zhang Get local grid boundaries
13060f0b76eSHong Zhang */
1319566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
13260f0b76eSHong Zhang
13360f0b76eSHong Zhang stencil[0].k = 0;
13460f0b76eSHong Zhang stencil[1].k = 0;
13560f0b76eSHong Zhang stencil[2].k = 0;
13660f0b76eSHong Zhang stencil[3].k = 0;
13760f0b76eSHong Zhang stencil[4].k = 0;
13860f0b76eSHong Zhang stencil[5].k = 0;
13960f0b76eSHong Zhang rowstencil.k = 0;
14060f0b76eSHong Zhang /*
14160f0b76eSHong Zhang Compute function over the locally owned part of the grid
14260f0b76eSHong Zhang */
14360f0b76eSHong Zhang for (j = ys; j < ys + ym; j++) {
14460f0b76eSHong Zhang stencil[0].j = j - 1;
14560f0b76eSHong Zhang stencil[1].j = j + 1;
14660f0b76eSHong Zhang stencil[2].j = j;
14760f0b76eSHong Zhang stencil[3].j = j;
14860f0b76eSHong Zhang stencil[4].j = j;
14960f0b76eSHong Zhang stencil[5].j = j;
1509371c9d4SSatish Balay rowstencil.k = 0;
1519371c9d4SSatish Balay rowstencil.j = j;
15260f0b76eSHong Zhang for (i = xs; i < xs + xm; i++) {
15360f0b76eSHong Zhang uc = u[j][i].u;
15460f0b76eSHong Zhang vc = u[j][i].v;
15560f0b76eSHong Zhang
15660f0b76eSHong Zhang /* uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
15760f0b76eSHong Zhang uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
15860f0b76eSHong Zhang
15960f0b76eSHong Zhang vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
16060f0b76eSHong Zhang vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
16160f0b76eSHong Zhang f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
16260f0b76eSHong Zhang
1639371c9d4SSatish Balay stencil[0].i = i;
1649371c9d4SSatish Balay stencil[0].c = 0;
1659371c9d4SSatish Balay entries[0] = appctx->D1 * sy;
1669371c9d4SSatish Balay stencil[1].i = i;
1679371c9d4SSatish Balay stencil[1].c = 0;
1689371c9d4SSatish Balay entries[1] = appctx->D1 * sy;
1699371c9d4SSatish Balay stencil[2].i = i - 1;
1709371c9d4SSatish Balay stencil[2].c = 0;
1719371c9d4SSatish Balay entries[2] = appctx->D1 * sx;
1729371c9d4SSatish Balay stencil[3].i = i + 1;
1739371c9d4SSatish Balay stencil[3].c = 0;
1749371c9d4SSatish Balay entries[3] = appctx->D1 * sx;
1759371c9d4SSatish Balay stencil[4].i = i;
1769371c9d4SSatish Balay stencil[4].c = 0;
1779371c9d4SSatish Balay entries[4] = -2.0 * appctx->D1 * (sx + sy) - vc * vc - appctx->gamma;
1789371c9d4SSatish Balay stencil[5].i = i;
1799371c9d4SSatish Balay stencil[5].c = 1;
1809371c9d4SSatish Balay entries[5] = -2.0 * uc * vc;
1819371c9d4SSatish Balay rowstencil.i = i;
1829371c9d4SSatish Balay rowstencil.c = 0;
18360f0b76eSHong Zhang
1849566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
18548a46eb9SPierre Jolivet if (appctx->aijpc) PetscCall(MatSetValuesStencil(BB, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
1869371c9d4SSatish Balay stencil[0].c = 1;
1879371c9d4SSatish Balay entries[0] = appctx->D2 * sy;
1889371c9d4SSatish Balay stencil[1].c = 1;
1899371c9d4SSatish Balay entries[1] = appctx->D2 * sy;
1909371c9d4SSatish Balay stencil[2].c = 1;
1919371c9d4SSatish Balay entries[2] = appctx->D2 * sx;
1929371c9d4SSatish Balay stencil[3].c = 1;
1939371c9d4SSatish Balay entries[3] = appctx->D2 * sx;
1949371c9d4SSatish Balay stencil[4].c = 1;
1959371c9d4SSatish Balay entries[4] = -2.0 * appctx->D2 * (sx + sy) + 2.0 * uc * vc - appctx->gamma - appctx->kappa;
1969371c9d4SSatish Balay stencil[5].c = 0;
1979371c9d4SSatish Balay entries[5] = vc * vc;
19860f0b76eSHong Zhang rowstencil.c = 1;
19960f0b76eSHong Zhang
2009566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
20148a46eb9SPierre Jolivet if (appctx->aijpc) PetscCall(MatSetValuesStencil(BB, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
20260f0b76eSHong Zhang /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
20360f0b76eSHong Zhang }
20460f0b76eSHong Zhang }
20560f0b76eSHong Zhang
20660f0b76eSHong Zhang /*
20760f0b76eSHong Zhang Restore vectors
20860f0b76eSHong Zhang */
2099566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(19.0 * xm * ym));
2109566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
2119566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
2129566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2139566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2149566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
21560f0b76eSHong Zhang if (appctx->aijpc) {
2169566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(BB, MAT_FINAL_ASSEMBLY));
2179566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(BB, MAT_FINAL_ASSEMBLY));
2189566063dSJacob Faibussowitsch PetscCall(MatSetOption(BB, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
21960f0b76eSHong Zhang }
2203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
22160f0b76eSHong Zhang }
22260f0b76eSHong Zhang
22360f0b76eSHong Zhang /*
22460f0b76eSHong Zhang IFunction - Evaluates implicit nonlinear function, xdot - F(x).
22560f0b76eSHong Zhang
22660f0b76eSHong Zhang Input Parameters:
22760f0b76eSHong Zhang . ts - the TS context
22860f0b76eSHong Zhang . U - input vector
22960f0b76eSHong Zhang . Udot - input vector
23060f0b76eSHong Zhang . ptr - optional user-defined context, as set by TSSetRHSFunction()
23160f0b76eSHong Zhang
23260f0b76eSHong Zhang Output Parameter:
23360f0b76eSHong Zhang . F - function vector
23460f0b76eSHong Zhang */
IFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void * ptr)235d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal ftime, Vec U, Vec Udot, Vec F, void *ptr)
236d71ae5a4SJacob Faibussowitsch {
23760f0b76eSHong Zhang AppCtx *appctx = (AppCtx *)ptr;
23860f0b76eSHong Zhang DM da;
23960f0b76eSHong Zhang PetscInt i, j, Mx, My, xs, ys, xm, ym;
24060f0b76eSHong Zhang PetscReal hx, hy, sx, sy;
24160f0b76eSHong Zhang PetscScalar uc, uxx, uyy, vc, vxx, vyy;
24260f0b76eSHong Zhang Field **u, **f, **udot;
24360f0b76eSHong Zhang Vec localU;
24460f0b76eSHong Zhang
24560f0b76eSHong Zhang PetscFunctionBegin;
2469566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
2479566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
2489566063dSJacob 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));
2499371c9d4SSatish Balay hx = 2.50 / (PetscReal)Mx;
2509371c9d4SSatish Balay sx = 1.0 / (hx * hx);
2519371c9d4SSatish Balay hy = 2.50 / (PetscReal)My;
2529371c9d4SSatish Balay sy = 1.0 / (hy * hy);
25360f0b76eSHong Zhang
25460f0b76eSHong Zhang /*
25560f0b76eSHong Zhang Scatter ghost points to local vector,using the 2-step process
25660f0b76eSHong Zhang DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
25760f0b76eSHong Zhang By placing code between these two statements, computations can be
25860f0b76eSHong Zhang done while messages are in transition.
25960f0b76eSHong Zhang */
2609566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
2619566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
26260f0b76eSHong Zhang
26360f0b76eSHong Zhang /*
26460f0b76eSHong Zhang Get pointers to vector data
26560f0b76eSHong Zhang */
2669566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u));
2679566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f));
2689566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Udot, &udot));
26960f0b76eSHong Zhang
27060f0b76eSHong Zhang /*
27160f0b76eSHong Zhang Get local grid boundaries
27260f0b76eSHong Zhang */
2739566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
27460f0b76eSHong Zhang
27560f0b76eSHong Zhang /*
27660f0b76eSHong Zhang Compute function over the locally owned part of the grid
27760f0b76eSHong Zhang */
27860f0b76eSHong Zhang for (j = ys; j < ys + ym; j++) {
27960f0b76eSHong Zhang for (i = xs; i < xs + xm; i++) {
28060f0b76eSHong Zhang uc = u[j][i].u;
28160f0b76eSHong Zhang uxx = (-2.0 * uc + u[j][i - 1].u + u[j][i + 1].u) * sx;
28260f0b76eSHong Zhang uyy = (-2.0 * uc + u[j - 1][i].u + u[j + 1][i].u) * sy;
28360f0b76eSHong Zhang vc = u[j][i].v;
28460f0b76eSHong Zhang vxx = (-2.0 * vc + u[j][i - 1].v + u[j][i + 1].v) * sx;
28560f0b76eSHong Zhang vyy = (-2.0 * vc + u[j - 1][i].v + u[j + 1][i].v) * sy;
28660f0b76eSHong Zhang f[j][i].u = udot[j][i].u - (appctx->D1 * (uxx + uyy) - uc * vc * vc + appctx->gamma * (1.0 - uc));
28760f0b76eSHong Zhang f[j][i].v = udot[j][i].v - (appctx->D2 * (vxx + vyy) + uc * vc * vc - (appctx->gamma + appctx->kappa) * vc);
28860f0b76eSHong Zhang }
28960f0b76eSHong Zhang }
2909566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym));
29160f0b76eSHong Zhang
29260f0b76eSHong Zhang /*
29360f0b76eSHong Zhang Restore vectors
29460f0b76eSHong Zhang */
2959566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
2969566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f));
2979566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Udot, &udot));
2989566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
2993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
30060f0b76eSHong Zhang }
30160f0b76eSHong Zhang
IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat BB,PetscCtx ctx)302*2a8381b2SBarry Smith PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat BB, PetscCtx ctx)
303d71ae5a4SJacob Faibussowitsch {
30460f0b76eSHong Zhang AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
30560f0b76eSHong Zhang DM da;
30660f0b76eSHong Zhang PetscInt i, j, Mx, My, xs, ys, xm, ym;
30760f0b76eSHong Zhang PetscReal hx, hy, sx, sy;
30860f0b76eSHong Zhang PetscScalar uc, vc;
30960f0b76eSHong Zhang Field **u;
31060f0b76eSHong Zhang Vec localU;
31160f0b76eSHong Zhang MatStencil stencil[6], rowstencil;
31260f0b76eSHong Zhang PetscScalar entries[6];
31360f0b76eSHong Zhang
31460f0b76eSHong Zhang PetscFunctionBegin;
3159566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
3169566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
3179566063dSJacob 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));
31860f0b76eSHong Zhang
3199371c9d4SSatish Balay hx = 2.50 / (PetscReal)Mx;
3209371c9d4SSatish Balay sx = 1.0 / (hx * hx);
3219371c9d4SSatish Balay hy = 2.50 / (PetscReal)My;
3229371c9d4SSatish Balay sy = 1.0 / (hy * hy);
32360f0b76eSHong Zhang
32460f0b76eSHong Zhang /*
32560f0b76eSHong Zhang Scatter ghost points to local vector,using the 2-step process
32660f0b76eSHong Zhang DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
32760f0b76eSHong Zhang By placing code between these two statements, computations can be
32860f0b76eSHong Zhang done while messages are in transition.
32960f0b76eSHong Zhang */
3309566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
3319566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
33260f0b76eSHong Zhang
33360f0b76eSHong Zhang /*
33460f0b76eSHong Zhang Get pointers to vector data
33560f0b76eSHong Zhang */
3369566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u));
33760f0b76eSHong Zhang
33860f0b76eSHong Zhang /*
33960f0b76eSHong Zhang Get local grid boundaries
34060f0b76eSHong Zhang */
3419566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
34260f0b76eSHong Zhang
34360f0b76eSHong Zhang stencil[0].k = 0;
34460f0b76eSHong Zhang stencil[1].k = 0;
34560f0b76eSHong Zhang stencil[2].k = 0;
34660f0b76eSHong Zhang stencil[3].k = 0;
34760f0b76eSHong Zhang stencil[4].k = 0;
34860f0b76eSHong Zhang stencil[5].k = 0;
34960f0b76eSHong Zhang rowstencil.k = 0;
35060f0b76eSHong Zhang /*
35160f0b76eSHong Zhang Compute function over the locally owned part of the grid
35260f0b76eSHong Zhang */
35360f0b76eSHong Zhang for (j = ys; j < ys + ym; j++) {
35460f0b76eSHong Zhang stencil[0].j = j - 1;
35560f0b76eSHong Zhang stencil[1].j = j + 1;
35660f0b76eSHong Zhang stencil[2].j = j;
35760f0b76eSHong Zhang stencil[3].j = j;
35860f0b76eSHong Zhang stencil[4].j = j;
35960f0b76eSHong Zhang stencil[5].j = j;
3609371c9d4SSatish Balay rowstencil.k = 0;
3619371c9d4SSatish Balay rowstencil.j = j;
36260f0b76eSHong Zhang for (i = xs; i < xs + xm; i++) {
36360f0b76eSHong Zhang uc = u[j][i].u;
36460f0b76eSHong Zhang vc = u[j][i].v;
36560f0b76eSHong Zhang
36660f0b76eSHong Zhang /* uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
36760f0b76eSHong Zhang uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
36860f0b76eSHong Zhang
36960f0b76eSHong Zhang vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
37060f0b76eSHong Zhang vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
37160f0b76eSHong Zhang f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
37260f0b76eSHong Zhang
3739371c9d4SSatish Balay stencil[0].i = i;
3749371c9d4SSatish Balay stencil[0].c = 0;
3759371c9d4SSatish Balay entries[0] = -appctx->D1 * sy;
3769371c9d4SSatish Balay stencil[1].i = i;
3779371c9d4SSatish Balay stencil[1].c = 0;
3789371c9d4SSatish Balay entries[1] = -appctx->D1 * sy;
3799371c9d4SSatish Balay stencil[2].i = i - 1;
3809371c9d4SSatish Balay stencil[2].c = 0;
3819371c9d4SSatish Balay entries[2] = -appctx->D1 * sx;
3829371c9d4SSatish Balay stencil[3].i = i + 1;
3839371c9d4SSatish Balay stencil[3].c = 0;
3849371c9d4SSatish Balay entries[3] = -appctx->D1 * sx;
3859371c9d4SSatish Balay stencil[4].i = i;
3869371c9d4SSatish Balay stencil[4].c = 0;
3879371c9d4SSatish Balay entries[4] = 2.0 * appctx->D1 * (sx + sy) + vc * vc + appctx->gamma + a;
3889371c9d4SSatish Balay stencil[5].i = i;
3899371c9d4SSatish Balay stencil[5].c = 1;
3909371c9d4SSatish Balay entries[5] = 2.0 * uc * vc;
3919371c9d4SSatish Balay rowstencil.i = i;
3929371c9d4SSatish Balay rowstencil.c = 0;
39360f0b76eSHong Zhang
3949566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
39548a46eb9SPierre Jolivet if (appctx->aijpc) PetscCall(MatSetValuesStencil(BB, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
3969371c9d4SSatish Balay stencil[0].c = 1;
3979371c9d4SSatish Balay entries[0] = -appctx->D2 * sy;
3989371c9d4SSatish Balay stencil[1].c = 1;
3999371c9d4SSatish Balay entries[1] = -appctx->D2 * sy;
4009371c9d4SSatish Balay stencil[2].c = 1;
4019371c9d4SSatish Balay entries[2] = -appctx->D2 * sx;
4029371c9d4SSatish Balay stencil[3].c = 1;
4039371c9d4SSatish Balay entries[3] = -appctx->D2 * sx;
4049371c9d4SSatish Balay stencil[4].c = 1;
4059371c9d4SSatish Balay entries[4] = 2.0 * appctx->D2 * (sx + sy) - 2.0 * uc * vc + appctx->gamma + appctx->kappa + a;
4069371c9d4SSatish Balay stencil[5].c = 0;
4079371c9d4SSatish Balay entries[5] = -vc * vc;
40860f0b76eSHong Zhang rowstencil.c = 1;
40960f0b76eSHong Zhang
4109566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
41148a46eb9SPierre Jolivet if (appctx->aijpc) PetscCall(MatSetValuesStencil(BB, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
41260f0b76eSHong Zhang /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
41360f0b76eSHong Zhang }
41460f0b76eSHong Zhang }
41560f0b76eSHong Zhang
41660f0b76eSHong Zhang /*
41760f0b76eSHong Zhang Restore vectors
41860f0b76eSHong Zhang */
4199566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(19.0 * xm * ym));
4209566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
4219566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
4229566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4239566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
4249566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
42560f0b76eSHong Zhang if (appctx->aijpc) {
4269566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(BB, MAT_FINAL_ASSEMBLY));
4279566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(BB, MAT_FINAL_ASSEMBLY));
4289566063dSJacob Faibussowitsch PetscCall(MatSetOption(BB, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
42960f0b76eSHong Zhang }
4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
43160f0b76eSHong Zhang }
432