xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex5adj_mf.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
160f0b76eSHong Zhang static char help[] = "Demonstrates adjoint sensitivity analysis for Reaction-Diffusion Equations.\n";
260f0b76eSHong Zhang 
360f0b76eSHong Zhang /*
460f0b76eSHong Zhang   See ex5.c for details on the equation.
560f0b76eSHong Zhang   This code demonestrates the TSAdjoint interface to a system of time-dependent partial differential equations.
660f0b76eSHong Zhang   It computes the sensitivity of a component in the final solution, which locates in the center of the 2D domain, w.r.t. the initial conditions.
760f0b76eSHong Zhang   The user does not need to provide any additional functions. The required functions in the original simulation are reused in the adjoint run.
860f0b76eSHong Zhang 
960f0b76eSHong Zhang   Runtime options:
1060f0b76eSHong Zhang     -forwardonly  - run the forward simulation without adjoint
1160f0b76eSHong Zhang     -implicitform - provide IFunction and IJacobian to TS, if not set, RHSFunction and RHSJacobian will be used
1260f0b76eSHong Zhang     -aijpc        - set the preconditioner matrix to be aij (the Jacobian matrix can be of a different type such as ELL)
1360f0b76eSHong Zhang */
1460f0b76eSHong Zhang 
1560f0b76eSHong Zhang #include "reaction_diffusion.h"
1660f0b76eSHong Zhang #include <petscdm.h>
1760f0b76eSHong Zhang #include <petscdmda.h>
1860f0b76eSHong Zhang 
1960f0b76eSHong Zhang /* Matrix free support */
2060f0b76eSHong Zhang typedef struct {
2160f0b76eSHong Zhang   PetscReal time;
2260f0b76eSHong Zhang   Vec       U;
2360f0b76eSHong Zhang   Vec       Udot;
2460f0b76eSHong Zhang   PetscReal shift;
2560f0b76eSHong Zhang   AppCtx   *appctx;
2660f0b76eSHong Zhang   TS        ts;
2760f0b76eSHong Zhang } MCtx;
2860f0b76eSHong Zhang 
2960f0b76eSHong Zhang /*
3060f0b76eSHong Zhang    User-defined routines
3160f0b76eSHong Zhang */
3260f0b76eSHong Zhang PetscErrorCode InitialConditions(DM, Vec);
3360f0b76eSHong Zhang PetscErrorCode RHSJacobianShell(TS, PetscReal, Vec, Mat, Mat, void *);
3460f0b76eSHong Zhang PetscErrorCode IJacobianShell(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
3560f0b76eSHong Zhang 
36d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeLambda(DM da, Vec lambda, PetscReal x, PetscReal y)
37d71ae5a4SJacob Faibussowitsch {
3860f0b76eSHong Zhang   PetscInt i, j, Mx, My, xs, ys, xm, ym;
3960f0b76eSHong Zhang   Field  **l;
4060f0b76eSHong Zhang   PetscFunctionBegin;
4160f0b76eSHong Zhang 
429566063dSJacob 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));
4360f0b76eSHong Zhang   /* locate the global i index for x and j index for y */
4460f0b76eSHong Zhang   i = (PetscInt)(x * (Mx - 1));
4560f0b76eSHong Zhang   j = (PetscInt)(y * (My - 1));
469566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
4760f0b76eSHong Zhang 
4860f0b76eSHong Zhang   if (xs <= i && i < xs + xm && ys <= j && j < ys + ym) {
4960f0b76eSHong Zhang     /* the i,j vertex is on this process */
509566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da, lambda, &l));
5160f0b76eSHong Zhang     l[j][i].u = 1.0;
5260f0b76eSHong Zhang     l[j][i].v = 1.0;
539566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da, lambda, &l));
5460f0b76eSHong Zhang   }
55*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5660f0b76eSHong Zhang }
5760f0b76eSHong Zhang 
58d71ae5a4SJacob Faibussowitsch static PetscErrorCode MyRHSMatMultTranspose(Mat A_shell, Vec X, Vec Y)
59d71ae5a4SJacob Faibussowitsch {
6060f0b76eSHong Zhang   MCtx       *mctx;
6160f0b76eSHong Zhang   AppCtx     *appctx;
6260f0b76eSHong Zhang   DM          da;
6360f0b76eSHong Zhang   PetscInt    i, j, Mx, My, xs, ys, xm, ym;
6460f0b76eSHong Zhang   PetscReal   hx, hy, sx, sy;
6560f0b76eSHong Zhang   PetscScalar uc, uxx, uyy, vc, vxx, vyy, ucb, vcb;
6660f0b76eSHong Zhang   Field     **u, **x, **y;
6760f0b76eSHong Zhang   Vec         localX;
6860f0b76eSHong Zhang 
6960f0b76eSHong Zhang   PetscFunctionBeginUser;
709566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A_shell, &mctx));
7160f0b76eSHong Zhang   appctx = mctx->appctx;
729566063dSJacob Faibussowitsch   PetscCall(TSGetDM(mctx->ts, &da));
739566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
749566063dSJacob 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));
759371c9d4SSatish Balay   hx = 2.50 / (PetscReal)Mx;
769371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
779371c9d4SSatish Balay   hy = 2.50 / (PetscReal)My;
789371c9d4SSatish Balay   sy = 1.0 / (hy * hy);
7960f0b76eSHong Zhang 
8060f0b76eSHong Zhang   /* Scatter ghost points to local vector,using the 2-step process
8160f0b76eSHong Zhang      DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
8260f0b76eSHong Zhang      By placing code between these two statements, computations can be
8360f0b76eSHong Zhang      done while messages are in transition. */
849566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
859566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
8660f0b76eSHong Zhang 
8760f0b76eSHong Zhang   /* Get pointers to vector data */
889566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localX, &x));
899566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, mctx->U, &u));
909566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Y, &y));
9160f0b76eSHong Zhang 
9260f0b76eSHong Zhang   /* Get local grid boundaries */
939566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
9460f0b76eSHong Zhang 
9560f0b76eSHong Zhang   /* Compute function over the locally owned part of the grid */
9660f0b76eSHong Zhang   for (j = ys; j < ys + ym; j++) {
9760f0b76eSHong Zhang     for (i = xs; i < xs + xm; i++) {
9860f0b76eSHong Zhang       uc        = u[j][i].u;
9960f0b76eSHong Zhang       ucb       = x[j][i].u;
10060f0b76eSHong Zhang       uxx       = (-2.0 * ucb + x[j][i - 1].u + x[j][i + 1].u) * sx;
10160f0b76eSHong Zhang       uyy       = (-2.0 * ucb + x[j - 1][i].u + x[j + 1][i].u) * sy;
10260f0b76eSHong Zhang       vc        = u[j][i].v;
10360f0b76eSHong Zhang       vcb       = x[j][i].v;
10460f0b76eSHong Zhang       vxx       = (-2.0 * vcb + x[j][i - 1].v + x[j][i + 1].v) * sx;
10560f0b76eSHong Zhang       vyy       = (-2.0 * vcb + x[j - 1][i].v + x[j + 1][i].v) * sy;
10660f0b76eSHong Zhang       y[j][i].u = appctx->D1 * (uxx + uyy) - ucb * (vc * vc + appctx->gamma) + vc * vc * vcb;
10760f0b76eSHong Zhang       y[j][i].v = appctx->D2 * (vxx + vyy) - 2.0 * uc * vc * ucb + (2.0 * uc * vc - appctx->gamma - appctx->kappa) * vcb;
10860f0b76eSHong Zhang     }
10960f0b76eSHong Zhang   }
1109566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localX, &x));
1119566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, mctx->U, &u));
1129566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Y, &y));
1139566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
114*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11560f0b76eSHong Zhang }
11660f0b76eSHong Zhang 
117d71ae5a4SJacob Faibussowitsch static PetscErrorCode MyIMatMultTranspose(Mat A_shell, Vec X, Vec Y)
118d71ae5a4SJacob Faibussowitsch {
11960f0b76eSHong Zhang   MCtx       *mctx;
12060f0b76eSHong Zhang   AppCtx     *appctx;
12160f0b76eSHong Zhang   DM          da;
12260f0b76eSHong Zhang   PetscInt    i, j, Mx, My, xs, ys, xm, ym;
12360f0b76eSHong Zhang   PetscReal   hx, hy, sx, sy;
12460f0b76eSHong Zhang   PetscScalar uc, uxx, uyy, vc, vxx, vyy, ucb, vcb;
12560f0b76eSHong Zhang   Field     **u, **x, **y;
12660f0b76eSHong Zhang   Vec         localX;
12760f0b76eSHong Zhang 
12860f0b76eSHong Zhang   PetscFunctionBeginUser;
1299566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A_shell, &mctx));
13060f0b76eSHong Zhang   appctx = mctx->appctx;
1319566063dSJacob Faibussowitsch   PetscCall(TSGetDM(mctx->ts, &da));
1329566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
1339566063dSJacob 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));
1349371c9d4SSatish Balay   hx = 2.50 / (PetscReal)Mx;
1359371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
1369371c9d4SSatish Balay   hy = 2.50 / (PetscReal)My;
1379371c9d4SSatish Balay   sy = 1.0 / (hy * hy);
13860f0b76eSHong Zhang 
13960f0b76eSHong Zhang   /* Scatter ghost points to local vector,using the 2-step process
14060f0b76eSHong Zhang      DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
14160f0b76eSHong Zhang      By placing code between these two statements, computations can be
14260f0b76eSHong Zhang      done while messages are in transition. */
1439566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
1449566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
14560f0b76eSHong Zhang 
14660f0b76eSHong Zhang   /* Get pointers to vector data */
1479566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localX, &x));
1489566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, mctx->U, &u));
1499566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Y, &y));
15060f0b76eSHong Zhang 
15160f0b76eSHong Zhang   /* Get local grid boundaries */
1529566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
15360f0b76eSHong Zhang 
15460f0b76eSHong Zhang   /* Compute function over the locally owned part of the grid */
15560f0b76eSHong Zhang   for (j = ys; j < ys + ym; j++) {
15660f0b76eSHong Zhang     for (i = xs; i < xs + xm; i++) {
15760f0b76eSHong Zhang       uc        = u[j][i].u;
15860f0b76eSHong Zhang       ucb       = x[j][i].u;
15960f0b76eSHong Zhang       uxx       = (-2.0 * ucb + x[j][i - 1].u + x[j][i + 1].u) * sx;
16060f0b76eSHong Zhang       uyy       = (-2.0 * ucb + x[j - 1][i].u + x[j + 1][i].u) * sy;
16160f0b76eSHong Zhang       vc        = u[j][i].v;
16260f0b76eSHong Zhang       vcb       = x[j][i].v;
16360f0b76eSHong Zhang       vxx       = (-2.0 * vcb + x[j][i - 1].v + x[j][i + 1].v) * sx;
16460f0b76eSHong Zhang       vyy       = (-2.0 * vcb + x[j - 1][i].v + x[j + 1][i].v) * sy;
16560f0b76eSHong Zhang       y[j][i].u = appctx->D1 * (uxx + uyy) - ucb * (vc * vc + appctx->gamma) + vc * vc * vcb;
16660f0b76eSHong Zhang       y[j][i].v = appctx->D2 * (vxx + vyy) - 2.0 * uc * vc * ucb + (2.0 * uc * vc - appctx->gamma - appctx->kappa) * vcb;
16760f0b76eSHong Zhang       y[j][i].u = mctx->shift * ucb - y[j][i].u;
16860f0b76eSHong Zhang       y[j][i].v = mctx->shift * vcb - y[j][i].v;
16960f0b76eSHong Zhang     }
17060f0b76eSHong Zhang   }
1719566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localX, &x));
1729566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, mctx->U, &u));
1739566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Y, &y));
1749566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
175*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17660f0b76eSHong Zhang }
17760f0b76eSHong Zhang 
178d71ae5a4SJacob Faibussowitsch static PetscErrorCode MyIMatMult(Mat A_shell, Vec X, Vec Y)
179d71ae5a4SJacob Faibussowitsch {
18060f0b76eSHong Zhang   MCtx       *mctx;
18160f0b76eSHong Zhang   AppCtx     *appctx;
18260f0b76eSHong Zhang   DM          da;
18360f0b76eSHong Zhang   PetscInt    i, j, Mx, My, xs, ys, xm, ym;
18460f0b76eSHong Zhang   PetscReal   hx, hy, sx, sy;
18560f0b76eSHong Zhang   PetscScalar uc, uxx, uyy, vc, vxx, vyy, ucb, vcb;
18660f0b76eSHong Zhang   Field     **u, **x, **y;
18760f0b76eSHong Zhang   Vec         localX;
18860f0b76eSHong Zhang 
18960f0b76eSHong Zhang   PetscFunctionBeginUser;
1909566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A_shell, &mctx));
19160f0b76eSHong Zhang   appctx = mctx->appctx;
1929566063dSJacob Faibussowitsch   PetscCall(TSGetDM(mctx->ts, &da));
1939566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localX));
1949566063dSJacob 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));
1959371c9d4SSatish Balay   hx = 2.50 / (PetscReal)Mx;
1969371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
1979371c9d4SSatish Balay   hy = 2.50 / (PetscReal)My;
1989371c9d4SSatish Balay   sy = 1.0 / (hy * hy);
19960f0b76eSHong Zhang 
20060f0b76eSHong Zhang   /* Scatter ghost points to local vector,using the 2-step process
20160f0b76eSHong Zhang      DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
20260f0b76eSHong Zhang      By placing code between these two statements, computations can be
20360f0b76eSHong Zhang      done while messages are in transition. */
2049566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, localX));
2059566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, localX));
20660f0b76eSHong Zhang 
20760f0b76eSHong Zhang   /* Get pointers to vector data */
2089566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localX, &x));
2099566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, mctx->U, &u));
2109566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, Y, &y));
21160f0b76eSHong Zhang 
21260f0b76eSHong Zhang   /* Get local grid boundaries */
2139566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
21460f0b76eSHong Zhang 
21560f0b76eSHong Zhang   /* Compute function over the locally owned part of the grid */
21660f0b76eSHong Zhang   for (j = ys; j < ys + ym; j++) {
21760f0b76eSHong Zhang     for (i = xs; i < xs + xm; i++) {
21860f0b76eSHong Zhang       uc        = u[j][i].u;
21960f0b76eSHong Zhang       ucb       = x[j][i].u;
22060f0b76eSHong Zhang       uxx       = (-2.0 * ucb + x[j][i - 1].u + x[j][i + 1].u) * sx;
22160f0b76eSHong Zhang       uyy       = (-2.0 * ucb + x[j - 1][i].u + x[j + 1][i].u) * sy;
22260f0b76eSHong Zhang       vc        = u[j][i].v;
22360f0b76eSHong Zhang       vcb       = x[j][i].v;
22460f0b76eSHong Zhang       vxx       = (-2.0 * vcb + x[j][i - 1].v + x[j][i + 1].v) * sx;
22560f0b76eSHong Zhang       vyy       = (-2.0 * vcb + x[j - 1][i].v + x[j + 1][i].v) * sy;
22660f0b76eSHong Zhang       y[j][i].u = appctx->D1 * (uxx + uyy) - (vc * vc + appctx->gamma) * ucb - 2.0 * uc * vc * vcb;
22760f0b76eSHong Zhang       y[j][i].v = appctx->D2 * (vxx + vyy) + vc * vc * ucb + (2.0 * uc * vc - appctx->gamma - appctx->kappa) * vcb;
22860f0b76eSHong Zhang       y[j][i].u = mctx->shift * ucb - y[j][i].u;
22960f0b76eSHong Zhang       y[j][i].v = mctx->shift * vcb - y[j][i].v;
23060f0b76eSHong Zhang     }
23160f0b76eSHong Zhang   }
2329566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localX, &x));
2339566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, mctx->U, &u));
2349566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, Y, &y));
2359566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localX));
236*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23760f0b76eSHong Zhang }
23860f0b76eSHong Zhang 
239d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
240d71ae5a4SJacob Faibussowitsch {
24160f0b76eSHong Zhang   TS             ts; /* ODE integrator */
24260f0b76eSHong Zhang   Vec            x;  /* solution */
24360f0b76eSHong Zhang   DM             da;
24460f0b76eSHong Zhang   AppCtx         appctx;
24560f0b76eSHong Zhang   MCtx           mctx;
24660f0b76eSHong Zhang   Vec            lambda[1];
24760f0b76eSHong Zhang   PetscBool      forwardonly = PETSC_FALSE, implicitform = PETSC_TRUE, mf = PETSC_FALSE;
24860f0b76eSHong Zhang   PetscLogDouble v1, v2;
24960f0b76eSHong Zhang #if defined(PETSC_USE_LOG)
25060f0b76eSHong Zhang   PetscLogStage stage;
25160f0b76eSHong Zhang #endif
25260f0b76eSHong Zhang 
253327415f7SBarry Smith   PetscFunctionBeginUser;
2549566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
2559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-forwardonly", &forwardonly, NULL));
2569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &implicitform, NULL));
25760f0b76eSHong Zhang   appctx.aijpc = PETSC_FALSE;
2589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-aijpc", &appctx.aijpc, NULL));
2599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-mf", &mf, NULL));
26060f0b76eSHong Zhang 
26160f0b76eSHong Zhang   appctx.D1    = 8.0e-5;
26260f0b76eSHong Zhang   appctx.D2    = 4.0e-5;
26360f0b76eSHong Zhang   appctx.gamma = .024;
26460f0b76eSHong Zhang   appctx.kappa = .06;
26560f0b76eSHong Zhang 
266*3ba16761SJacob Faibussowitsch   PetscCall(PetscLogStageRegister("MyAdjoint", &stage));
26760f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26860f0b76eSHong Zhang      Create distributed array (DMDA) to manage parallel grid and vectors
26960f0b76eSHong Zhang   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2709566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_PERIODIC, DMDA_STENCIL_STAR, 64, 64, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, NULL, &da));
2719566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
2729566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
2739566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "u"));
2749566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 1, "v"));
27560f0b76eSHong Zhang 
27660f0b76eSHong Zhang   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
27760f0b76eSHong Zhang      Extract global vectors from DMDA; then duplicate for remaining
27860f0b76eSHong Zhang      vectors that are the same types
27960f0b76eSHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2809566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
28160f0b76eSHong Zhang 
28260f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
28360f0b76eSHong Zhang      Create timestepping solver context
28460f0b76eSHong Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2859566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2869566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
2879566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
2889566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
28960f0b76eSHong Zhang   if (!implicitform) {
2909566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts, TSRK));
2919566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &appctx));
2929566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobian, &appctx));
29360f0b76eSHong Zhang   } else {
2949566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts, TSCN));
2959566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction, &appctx));
29660f0b76eSHong Zhang     if (appctx.aijpc) {
29760f0b76eSHong Zhang       Mat A, B;
29860f0b76eSHong Zhang 
2999566063dSJacob Faibussowitsch       PetscCall(DMSetMatType(da, MATSELL));
3009566063dSJacob Faibussowitsch       PetscCall(DMCreateMatrix(da, &A));
3019566063dSJacob Faibussowitsch       PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B));
30260f0b76eSHong Zhang       /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */
3039566063dSJacob Faibussowitsch       PetscCall(TSSetIJacobian(ts, A, B, IJacobian, &appctx));
3049566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A));
3059566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B));
30660f0b76eSHong Zhang     } else {
3079566063dSJacob Faibussowitsch       PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobian, &appctx));
30860f0b76eSHong Zhang     }
30960f0b76eSHong Zhang   }
31060f0b76eSHong Zhang 
31160f0b76eSHong Zhang   if (mf) {
31260f0b76eSHong Zhang     PetscInt xm, ym, Mx, My, dof;
31360f0b76eSHong Zhang     mctx.ts     = ts;
31460f0b76eSHong Zhang     mctx.appctx = &appctx;
3159566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &mctx.U));
31660f0b76eSHong Zhang     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31760f0b76eSHong Zhang       Create matrix free context
31860f0b76eSHong Zhang       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3199566063dSJacob Faibussowitsch     PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, &dof, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
3209566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(da, NULL, NULL, NULL, &xm, &ym, NULL));
3219566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PETSC_COMM_WORLD, dof * xm * ym, PETSC_DETERMINE, dof * Mx * My, dof * Mx * My, &mctx, &appctx.A));
3229566063dSJacob Faibussowitsch     PetscCall(MatShellSetOperation(appctx.A, MATOP_MULT_TRANSPOSE, (void (*)(void))MyRHSMatMultTranspose));
32360f0b76eSHong Zhang     if (!implicitform) { /* for explicit methods only */
3249566063dSJacob Faibussowitsch       PetscCall(TSSetRHSJacobian(ts, appctx.A, appctx.A, RHSJacobianShell, &appctx));
32560f0b76eSHong Zhang     } else {
3269566063dSJacob Faibussowitsch       /* PetscCall(VecDuplicate(appctx.U,&mctx.Udot)); */
3279566063dSJacob Faibussowitsch       PetscCall(MatShellSetOperation(appctx.A, MATOP_MULT, (void (*)(void))MyIMatMult));
3289566063dSJacob Faibussowitsch       PetscCall(MatShellSetOperation(appctx.A, MATOP_MULT_TRANSPOSE, (void (*)(void))MyIMatMultTranspose));
3299566063dSJacob Faibussowitsch       PetscCall(TSSetIJacobian(ts, appctx.A, appctx.A, IJacobianShell, &appctx));
33060f0b76eSHong Zhang     }
33160f0b76eSHong Zhang   }
33260f0b76eSHong Zhang 
33360f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33460f0b76eSHong Zhang      Set initial conditions
33560f0b76eSHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3369566063dSJacob Faibussowitsch   PetscCall(InitialConditions(da, x));
3379566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, x));
33860f0b76eSHong Zhang 
33960f0b76eSHong Zhang   /*
34060f0b76eSHong Zhang     Have the TS save its trajectory so that TSAdjointSolve() may be used
34160f0b76eSHong Zhang   */
3429566063dSJacob Faibussowitsch   if (!forwardonly) PetscCall(TSSetSaveTrajectory(ts));
34360f0b76eSHong Zhang 
34460f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34560f0b76eSHong Zhang      Set solver options
34660f0b76eSHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3479566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 200.0));
3489566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.5));
3499566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
3509566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
35160f0b76eSHong Zhang 
3529566063dSJacob Faibussowitsch   PetscCall(PetscTime(&v1));
35360f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35460f0b76eSHong Zhang      Solve ODE system
35560f0b76eSHong Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3569566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
35760f0b76eSHong Zhang   if (!forwardonly) {
35860f0b76eSHong Zhang     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35960f0b76eSHong Zhang        Start the Adjoint model
36060f0b76eSHong Zhang        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3619566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x, &lambda[0]));
36260f0b76eSHong Zhang     /*   Reset initial conditions for the adjoint integration */
3639566063dSJacob Faibussowitsch     PetscCall(InitializeLambda(da, lambda[0], 0.5, 0.5));
3649566063dSJacob Faibussowitsch     PetscCall(TSSetCostGradients(ts, 1, lambda, NULL));
365*3ba16761SJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage));
3669566063dSJacob Faibussowitsch     PetscCall(TSAdjointSolve(ts));
367*3ba16761SJacob Faibussowitsch     PetscCall(PetscLogStagePop());
3689566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lambda[0]));
36960f0b76eSHong Zhang   }
3709566063dSJacob Faibussowitsch   PetscCall(PetscTime(&v2));
3719566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, " %.3lf ", v2 - v1));
37260f0b76eSHong Zhang 
37360f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37460f0b76eSHong Zhang      Free work space.  All PETSc objects should be destroyed when they
37560f0b76eSHong Zhang      are no longer needed.
37660f0b76eSHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3789566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3799566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
38060f0b76eSHong Zhang   if (mf) {
3819566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mctx.U));
3829566063dSJacob Faibussowitsch     /* PetscCall(VecDestroy(&mctx.Udot));*/
3839566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&appctx.A));
38460f0b76eSHong Zhang   }
3859566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
386b122ec5aSJacob Faibussowitsch   return 0;
38760f0b76eSHong Zhang }
38860f0b76eSHong Zhang 
38960f0b76eSHong Zhang /* ------------------------------------------------------------------- */
390d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobianShell(TS ts, PetscReal t, Vec U, Mat A, Mat BB, void *ctx)
391d71ae5a4SJacob Faibussowitsch {
39260f0b76eSHong Zhang   MCtx *mctx;
39360f0b76eSHong Zhang 
39460f0b76eSHong Zhang   PetscFunctionBegin;
3959566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &mctx));
3969566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, mctx->U));
397*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39860f0b76eSHong Zhang }
39960f0b76eSHong Zhang 
400d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobianShell(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat BB, void *ctx)
401d71ae5a4SJacob Faibussowitsch {
40260f0b76eSHong Zhang   MCtx *mctx;
40360f0b76eSHong Zhang 
40460f0b76eSHong Zhang   PetscFunctionBegin;
4059566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &mctx));
4069566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, mctx->U));
4079566063dSJacob Faibussowitsch   /* PetscCall(VecCopy(Udot,mctx->Udot)); */
40860f0b76eSHong Zhang   mctx->shift = a;
409*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41060f0b76eSHong Zhang }
41160f0b76eSHong Zhang 
41260f0b76eSHong Zhang /* ------------------------------------------------------------------- */
413d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(DM da, Vec U)
414d71ae5a4SJacob Faibussowitsch {
41560f0b76eSHong Zhang   PetscInt  i, j, xs, ys, xm, ym, Mx, My;
41660f0b76eSHong Zhang   Field   **u;
41760f0b76eSHong Zhang   PetscReal hx, hy, x, y;
41860f0b76eSHong Zhang 
41960f0b76eSHong Zhang   PetscFunctionBegin;
4209566063dSJacob 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));
42160f0b76eSHong Zhang 
42260f0b76eSHong Zhang   hx = 2.5 / (PetscReal)Mx;
42360f0b76eSHong Zhang   hy = 2.5 / (PetscReal)My;
42460f0b76eSHong Zhang 
42560f0b76eSHong Zhang   /*
42660f0b76eSHong Zhang      Get pointers to vector data
42760f0b76eSHong Zhang   */
4289566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, U, &u));
42960f0b76eSHong Zhang 
43060f0b76eSHong Zhang   /*
43160f0b76eSHong Zhang      Get local grid boundaries
43260f0b76eSHong Zhang   */
4339566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
43460f0b76eSHong Zhang 
43560f0b76eSHong Zhang   /*
43660f0b76eSHong Zhang      Compute function over the locally owned part of the grid
43760f0b76eSHong Zhang   */
43860f0b76eSHong Zhang   for (j = ys; j < ys + ym; j++) {
43960f0b76eSHong Zhang     y = j * hy;
44060f0b76eSHong Zhang     for (i = xs; i < xs + xm; i++) {
44160f0b76eSHong Zhang       x = i * hx;
44260f0b76eSHong Zhang       if ((1.0 <= x) && (x <= 1.5) && (1.0 <= y) && (y <= 1.5)) u[j][i].v = .25 * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * x), 2.0) * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * y), 2.0);
44360f0b76eSHong Zhang       else u[j][i].v = 0.0;
44460f0b76eSHong Zhang 
44560f0b76eSHong Zhang       u[j][i].u = 1.0 - 2.0 * u[j][i].v;
44660f0b76eSHong Zhang     }
44760f0b76eSHong Zhang   }
44860f0b76eSHong Zhang 
44960f0b76eSHong Zhang   /*
45060f0b76eSHong Zhang      Restore vectors
45160f0b76eSHong Zhang   */
4529566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, U, &u));
453*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45460f0b76eSHong Zhang }
45560f0b76eSHong Zhang 
45660f0b76eSHong Zhang /*TEST
45760f0b76eSHong Zhang 
45860f0b76eSHong Zhang    build:
45960f0b76eSHong Zhang       depends: reaction_diffusion.c
46060f0b76eSHong Zhang       requires: !complex !single
46160f0b76eSHong Zhang 
46260f0b76eSHong Zhang    test:
46360f0b76eSHong Zhang       requires: cams
464533251d3SHong Zhang       args: -mf -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -ts_trajectory_type memory -ts_trajectory_solution_only 0 -ts_trajectory_max_units_ram 6 -ts_trajectory_memory_type cams
46560f0b76eSHong Zhang TEST*/
466