xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex5adj_mf.c (revision 533251d39df67bc1a9e2b8f423bcf8f28505e286)
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 
3660f0b76eSHong Zhang PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y)
3760f0b76eSHong Zhang {
3860f0b76eSHong Zhang    PetscInt i,j,Mx,My,xs,ys,xm,ym;
3960f0b76eSHong Zhang    PetscErrorCode ierr;
4060f0b76eSHong Zhang    Field **l;
4160f0b76eSHong Zhang    PetscFunctionBegin;
4260f0b76eSHong Zhang 
4360f0b76eSHong Zhang    ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
4460f0b76eSHong Zhang    /* locate the global i index for x and j index for y */
4560f0b76eSHong Zhang    i = (PetscInt)(x*(Mx-1));
4660f0b76eSHong Zhang    j = (PetscInt)(y*(My-1));
4760f0b76eSHong Zhang    ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
4860f0b76eSHong Zhang 
4960f0b76eSHong Zhang    if (xs <= i && i < xs+xm && ys <= j && j < ys+ym) {
5060f0b76eSHong Zhang      /* the i,j vertex is on this process */
5160f0b76eSHong Zhang      ierr = DMDAVecGetArray(da,lambda,&l);CHKERRQ(ierr);
5260f0b76eSHong Zhang      l[j][i].u = 1.0;
5360f0b76eSHong Zhang      l[j][i].v = 1.0;
5460f0b76eSHong Zhang      ierr = DMDAVecRestoreArray(da,lambda,&l);CHKERRQ(ierr);
5560f0b76eSHong Zhang    }
5660f0b76eSHong Zhang    PetscFunctionReturn(0);
5760f0b76eSHong Zhang }
5860f0b76eSHong Zhang 
5960f0b76eSHong Zhang static PetscErrorCode MyRHSMatMultTranspose(Mat A_shell,Vec X,Vec Y)
6060f0b76eSHong Zhang {
6160f0b76eSHong Zhang   MCtx           *mctx;
6260f0b76eSHong Zhang   AppCtx         *appctx;
6360f0b76eSHong Zhang   DM             da;
6460f0b76eSHong Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
6560f0b76eSHong Zhang   PetscReal      hx,hy,sx,sy;
6660f0b76eSHong Zhang   PetscScalar    uc,uxx,uyy,vc,vxx,vyy,ucb,vcb;
6760f0b76eSHong Zhang   Field          **u,**x,**y;
6860f0b76eSHong Zhang   Vec            localX;
6960f0b76eSHong Zhang   PetscErrorCode ierr;
7060f0b76eSHong Zhang 
7160f0b76eSHong Zhang   PetscFunctionBeginUser;
7260f0b76eSHong Zhang   ierr = MatShellGetContext(A_shell,&mctx);CHKERRQ(ierr);
7360f0b76eSHong Zhang   appctx = mctx->appctx;
7460f0b76eSHong Zhang   ierr = TSGetDM(mctx->ts,&da);CHKERRQ(ierr);
7560f0b76eSHong Zhang   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
7660f0b76eSHong Zhang   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
7760f0b76eSHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
7860f0b76eSHong Zhang   hy = 2.50/(PetscReal)My; 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. */
8460f0b76eSHong Zhang   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
8560f0b76eSHong Zhang   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
8660f0b76eSHong Zhang 
8760f0b76eSHong Zhang   /* Get pointers to vector data */
8860f0b76eSHong Zhang   ierr = DMDAVecGetArrayRead(da,localX,&x);CHKERRQ(ierr);
8960f0b76eSHong Zhang   ierr = DMDAVecGetArrayRead(da,mctx->U,&u);CHKERRQ(ierr);
9060f0b76eSHong Zhang   ierr = DMDAVecGetArray(da,Y,&y);CHKERRQ(ierr);
9160f0b76eSHong Zhang 
9260f0b76eSHong Zhang   /* Get local grid boundaries */
9360f0b76eSHong Zhang   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
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   }
11060f0b76eSHong Zhang   ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr);
11160f0b76eSHong Zhang   ierr = DMDAVecRestoreArrayRead(da,mctx->U,&u);CHKERRQ(ierr);
11260f0b76eSHong Zhang   ierr = DMDAVecRestoreArray(da,Y,&y);CHKERRQ(ierr);
11360f0b76eSHong Zhang   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
11460f0b76eSHong Zhang   PetscFunctionReturn(0);
11560f0b76eSHong Zhang }
11660f0b76eSHong Zhang 
11760f0b76eSHong Zhang static PetscErrorCode MyIMatMultTranspose(Mat A_shell,Vec X,Vec Y)
11860f0b76eSHong Zhang {
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   PetscErrorCode ierr;
12860f0b76eSHong Zhang 
12960f0b76eSHong Zhang   PetscFunctionBeginUser;
13060f0b76eSHong Zhang   ierr = MatShellGetContext(A_shell,&mctx);CHKERRQ(ierr);
13160f0b76eSHong Zhang   appctx = mctx->appctx;
13260f0b76eSHong Zhang   ierr = TSGetDM(mctx->ts,&da);CHKERRQ(ierr);
13360f0b76eSHong Zhang   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
13460f0b76eSHong Zhang   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
13560f0b76eSHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
13660f0b76eSHong Zhang   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
13760f0b76eSHong Zhang 
13860f0b76eSHong Zhang   /* Scatter ghost points to local vector,using the 2-step process
13960f0b76eSHong Zhang      DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
14060f0b76eSHong Zhang      By placing code between these two statements, computations can be
14160f0b76eSHong Zhang      done while messages are in transition. */
14260f0b76eSHong Zhang   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
14360f0b76eSHong Zhang   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
14460f0b76eSHong Zhang 
14560f0b76eSHong Zhang   /* Get pointers to vector data */
14660f0b76eSHong Zhang   ierr = DMDAVecGetArrayRead(da,localX,&x);CHKERRQ(ierr);
14760f0b76eSHong Zhang   ierr = DMDAVecGetArrayRead(da,mctx->U,&u);CHKERRQ(ierr);
14860f0b76eSHong Zhang   ierr = DMDAVecGetArray(da,Y,&y);CHKERRQ(ierr);
14960f0b76eSHong Zhang 
15060f0b76eSHong Zhang   /* Get local grid boundaries */
15160f0b76eSHong Zhang   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
15260f0b76eSHong Zhang 
15360f0b76eSHong Zhang   /* Compute function over the locally owned part of the grid */
15460f0b76eSHong Zhang   for (j=ys; j<ys+ym; j++) {
15560f0b76eSHong Zhang     for (i=xs; i<xs+xm; i++) {
15660f0b76eSHong Zhang       uc        = u[j][i].u;
15760f0b76eSHong Zhang       ucb       = x[j][i].u;
15860f0b76eSHong Zhang       uxx       = (-2.0*ucb + x[j][i-1].u + x[j][i+1].u)*sx;
15960f0b76eSHong Zhang       uyy       = (-2.0*ucb + x[j-1][i].u + x[j+1][i].u)*sy;
16060f0b76eSHong Zhang       vc        = u[j][i].v;
16160f0b76eSHong Zhang       vcb       = x[j][i].v;
16260f0b76eSHong Zhang       vxx       = (-2.0*vcb + x[j][i-1].v + x[j][i+1].v)*sx;
16360f0b76eSHong Zhang       vyy       = (-2.0*vcb + x[j-1][i].v + x[j+1][i].v)*sy;
16460f0b76eSHong Zhang       y[j][i].u = appctx->D1*(uxx + uyy) - ucb*(vc*vc+appctx->gamma) + vc*vc*vcb;
16560f0b76eSHong Zhang       y[j][i].v = appctx->D2*(vxx + vyy) - 2.0*uc*vc*ucb + (2.0*uc*vc-appctx->gamma - appctx->kappa)*vcb;
16660f0b76eSHong Zhang       y[j][i].u = mctx->shift*ucb - y[j][i].u;
16760f0b76eSHong Zhang       y[j][i].v = mctx->shift*vcb - y[j][i].v;
16860f0b76eSHong Zhang     }
16960f0b76eSHong Zhang   }
17060f0b76eSHong Zhang   ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr);
17160f0b76eSHong Zhang   ierr = DMDAVecRestoreArrayRead(da,mctx->U,&u);CHKERRQ(ierr);
17260f0b76eSHong Zhang   ierr = DMDAVecRestoreArray(da,Y,&y);CHKERRQ(ierr);
17360f0b76eSHong Zhang   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
17460f0b76eSHong Zhang   PetscFunctionReturn(0);
17560f0b76eSHong Zhang }
17660f0b76eSHong Zhang 
17760f0b76eSHong Zhang static PetscErrorCode MyIMatMult(Mat A_shell,Vec X,Vec Y)
17860f0b76eSHong Zhang {
17960f0b76eSHong Zhang   MCtx           *mctx;
18060f0b76eSHong Zhang   AppCtx         *appctx;
18160f0b76eSHong Zhang   DM             da;
18260f0b76eSHong Zhang   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
18360f0b76eSHong Zhang   PetscReal      hx,hy,sx,sy;
18460f0b76eSHong Zhang   PetscScalar    uc,uxx,uyy,vc,vxx,vyy,ucb,vcb;
18560f0b76eSHong Zhang   Field          **u,**x,**y;
18660f0b76eSHong Zhang   Vec            localX;
18760f0b76eSHong Zhang   PetscErrorCode ierr;
18860f0b76eSHong Zhang 
18960f0b76eSHong Zhang   PetscFunctionBeginUser;
19060f0b76eSHong Zhang   ierr = MatShellGetContext(A_shell,&mctx);CHKERRQ(ierr);
19160f0b76eSHong Zhang   appctx = mctx->appctx;
19260f0b76eSHong Zhang   ierr = TSGetDM(mctx->ts,&da);CHKERRQ(ierr);
19360f0b76eSHong Zhang   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
19460f0b76eSHong Zhang   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
19560f0b76eSHong Zhang   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
19660f0b76eSHong Zhang   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
19760f0b76eSHong Zhang 
19860f0b76eSHong Zhang   /* Scatter ghost points to local vector,using the 2-step process
19960f0b76eSHong Zhang      DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
20060f0b76eSHong Zhang      By placing code between these two statements, computations can be
20160f0b76eSHong Zhang      done while messages are in transition. */
20260f0b76eSHong Zhang   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
20360f0b76eSHong Zhang   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
20460f0b76eSHong Zhang 
20560f0b76eSHong Zhang   /* Get pointers to vector data */
20660f0b76eSHong Zhang   ierr = DMDAVecGetArrayRead(da,localX,&x);CHKERRQ(ierr);
20760f0b76eSHong Zhang   ierr = DMDAVecGetArrayRead(da,mctx->U,&u);CHKERRQ(ierr);
20860f0b76eSHong Zhang   ierr = DMDAVecGetArray(da,Y,&y);CHKERRQ(ierr);
20960f0b76eSHong Zhang 
21060f0b76eSHong Zhang   /* Get local grid boundaries */
21160f0b76eSHong Zhang   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
21260f0b76eSHong Zhang 
21360f0b76eSHong Zhang   /* Compute function over the locally owned part of the grid */
21460f0b76eSHong Zhang   for (j=ys; j<ys+ym; j++) {
21560f0b76eSHong Zhang     for (i=xs; i<xs+xm; i++) {
21660f0b76eSHong Zhang       uc        = u[j][i].u;
21760f0b76eSHong Zhang       ucb       = x[j][i].u;
21860f0b76eSHong Zhang       uxx       = (-2.0*ucb + x[j][i-1].u + x[j][i+1].u)*sx;
21960f0b76eSHong Zhang       uyy       = (-2.0*ucb + x[j-1][i].u + x[j+1][i].u)*sy;
22060f0b76eSHong Zhang       vc        = u[j][i].v;
22160f0b76eSHong Zhang       vcb       = x[j][i].v;
22260f0b76eSHong Zhang       vxx       = (-2.0*vcb + x[j][i-1].v + x[j][i+1].v)*sx;
22360f0b76eSHong Zhang       vyy       = (-2.0*vcb + x[j-1][i].v + x[j+1][i].v)*sy;
22460f0b76eSHong Zhang       y[j][i].u = appctx->D1*(uxx + uyy) - (vc*vc+appctx->gamma)*ucb - 2.0*uc*vc*vcb;
22560f0b76eSHong Zhang       y[j][i].v = appctx->D2*(vxx + vyy) + vc*vc*ucb + (2.0*uc*vc-appctx->gamma-appctx->kappa)*vcb;
22660f0b76eSHong Zhang       y[j][i].u = mctx->shift*ucb - y[j][i].u;
22760f0b76eSHong Zhang       y[j][i].v = mctx->shift*vcb - y[j][i].v;
22860f0b76eSHong Zhang     }
22960f0b76eSHong Zhang   }
23060f0b76eSHong Zhang   ierr = DMDAVecRestoreArrayRead(da,localX,&x);CHKERRQ(ierr);
23160f0b76eSHong Zhang   ierr = DMDAVecRestoreArrayRead(da,mctx->U,&u);CHKERRQ(ierr);
23260f0b76eSHong Zhang   ierr = DMDAVecRestoreArray(da,Y,&y);CHKERRQ(ierr);
23360f0b76eSHong Zhang   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
23460f0b76eSHong Zhang   PetscFunctionReturn(0);
23560f0b76eSHong Zhang }
23660f0b76eSHong Zhang 
23760f0b76eSHong Zhang int main(int argc,char **argv)
23860f0b76eSHong Zhang {
23960f0b76eSHong Zhang   TS             ts;                  /* ODE integrator */
24060f0b76eSHong Zhang   Vec            x;                   /* solution */
24160f0b76eSHong Zhang   PetscErrorCode ierr;
24260f0b76eSHong Zhang   DM             da;
24360f0b76eSHong Zhang   AppCtx         appctx;
24460f0b76eSHong Zhang   MCtx           mctx;
24560f0b76eSHong Zhang   Vec            lambda[1];
24660f0b76eSHong Zhang   PetscBool      forwardonly=PETSC_FALSE,implicitform=PETSC_TRUE,mf = PETSC_FALSE;
24760f0b76eSHong Zhang   PetscLogDouble v1,v2;
24860f0b76eSHong Zhang #if defined(PETSC_USE_LOG)
24960f0b76eSHong Zhang   PetscLogStage  stage;
25060f0b76eSHong Zhang #endif
25160f0b76eSHong Zhang 
25260f0b76eSHong Zhang   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
25360f0b76eSHong Zhang   ierr = PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL);CHKERRQ(ierr);
25460f0b76eSHong Zhang   ierr = PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL);CHKERRQ(ierr);
25560f0b76eSHong Zhang   appctx.aijpc = PETSC_FALSE;
25660f0b76eSHong Zhang   ierr = PetscOptionsGetBool(NULL,NULL,"-aijpc",&appctx.aijpc,NULL);CHKERRQ(ierr);
25760f0b76eSHong Zhang   ierr = PetscOptionsGetBool(NULL,NULL,"-mf",&mf,NULL);CHKERRQ(ierr);
25860f0b76eSHong Zhang 
25960f0b76eSHong Zhang   appctx.D1    = 8.0e-5;
26060f0b76eSHong Zhang   appctx.D2    = 4.0e-5;
26160f0b76eSHong Zhang   appctx.gamma = .024;
26260f0b76eSHong Zhang   appctx.kappa = .06;
26360f0b76eSHong Zhang 
26460f0b76eSHong Zhang   PetscLogStageRegister("MyAdjoint", &stage);
26560f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26660f0b76eSHong Zhang      Create distributed array (DMDA) to manage parallel grid and vectors
26760f0b76eSHong Zhang   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
26860f0b76eSHong Zhang   ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,64,64,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da);CHKERRQ(ierr);
26960f0b76eSHong Zhang   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
27060f0b76eSHong Zhang   ierr = DMSetUp(da);CHKERRQ(ierr);
27160f0b76eSHong Zhang   ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr);
27260f0b76eSHong Zhang   ierr = DMDASetFieldName(da,1,"v");CHKERRQ(ierr);
27360f0b76eSHong Zhang 
27460f0b76eSHong Zhang   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
27560f0b76eSHong Zhang      Extract global vectors from DMDA; then duplicate for remaining
27660f0b76eSHong Zhang      vectors that are the same types
27760f0b76eSHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
27860f0b76eSHong Zhang   ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr);
27960f0b76eSHong Zhang 
28060f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
28160f0b76eSHong Zhang      Create timestepping solver context
28260f0b76eSHong Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
28360f0b76eSHong Zhang   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
28460f0b76eSHong Zhang   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
28560f0b76eSHong Zhang   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
28660f0b76eSHong Zhang   ierr = TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
28760f0b76eSHong Zhang   if (!implicitform) {
28860f0b76eSHong Zhang     ierr = TSSetType(ts,TSRK);CHKERRQ(ierr);
28960f0b76eSHong Zhang     ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&appctx);CHKERRQ(ierr);
29060f0b76eSHong Zhang     ierr = TSSetRHSJacobian(ts,NULL,NULL,RHSJacobian,&appctx);CHKERRQ(ierr);
29160f0b76eSHong Zhang   } else {
29260f0b76eSHong Zhang     ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
29360f0b76eSHong Zhang     ierr = TSSetIFunction(ts,NULL,IFunction,&appctx);CHKERRQ(ierr);
29460f0b76eSHong Zhang     if (appctx.aijpc) {
29560f0b76eSHong Zhang       Mat                    A,B;
29660f0b76eSHong Zhang 
29760f0b76eSHong Zhang       ierr = DMSetMatType(da,MATSELL);CHKERRQ(ierr);
29860f0b76eSHong Zhang       ierr = DMCreateMatrix(da,&A);CHKERRQ(ierr);
29960f0b76eSHong Zhang       ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
30060f0b76eSHong Zhang       /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */
30160f0b76eSHong Zhang       ierr = TSSetIJacobian(ts,A,B,IJacobian,&appctx);CHKERRQ(ierr);
30260f0b76eSHong Zhang       ierr = MatDestroy(&A);CHKERRQ(ierr);
30360f0b76eSHong Zhang       ierr = MatDestroy(&B);CHKERRQ(ierr);
30460f0b76eSHong Zhang     } else {
30560f0b76eSHong Zhang       ierr = TSSetIJacobian(ts,NULL,NULL,IJacobian,&appctx);CHKERRQ(ierr);
30660f0b76eSHong Zhang     }
30760f0b76eSHong Zhang   }
30860f0b76eSHong Zhang 
30960f0b76eSHong Zhang   if (mf) {
31060f0b76eSHong Zhang     PetscInt xm,ym,Mx,My,dof;
31160f0b76eSHong Zhang     mctx.ts = ts;
31260f0b76eSHong Zhang     mctx.appctx = &appctx;
31360f0b76eSHong Zhang     ierr = VecDuplicate(x,&mctx.U);CHKERRQ(ierr);
31460f0b76eSHong Zhang     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31560f0b76eSHong Zhang       Create matrix free context
31660f0b76eSHong Zhang       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
31760f0b76eSHong Zhang     ierr = 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);CHKERRQ(ierr);
31860f0b76eSHong Zhang     ierr = DMDAGetCorners(da,NULL,NULL,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
31960f0b76eSHong Zhang     ierr = MatCreateShell(PETSC_COMM_WORLD,dof*xm*ym,PETSC_DETERMINE,dof*Mx*My,dof*Mx*My,&mctx,&appctx.A);CHKERRQ(ierr);
32060f0b76eSHong Zhang     ierr = MatShellSetOperation(appctx.A,MATOP_MULT_TRANSPOSE,(void (*)(void))MyRHSMatMultTranspose);CHKERRQ(ierr);
32160f0b76eSHong Zhang     if (!implicitform) { /* for explicit methods only */
32260f0b76eSHong Zhang       ierr = TSSetRHSJacobian(ts,appctx.A,appctx.A,RHSJacobianShell,&appctx);CHKERRQ(ierr);
32360f0b76eSHong Zhang     } else {
32460f0b76eSHong Zhang       /* ierr = VecDuplicate(appctx.U,&mctx.Udot);CHKERRQ(ierr); */
32560f0b76eSHong Zhang       ierr = MatShellSetOperation(appctx.A,MATOP_MULT,(void (*)(void))MyIMatMult);CHKERRQ(ierr);
32660f0b76eSHong Zhang       ierr = MatShellSetOperation(appctx.A,MATOP_MULT_TRANSPOSE,(void (*)(void))MyIMatMultTranspose);CHKERRQ(ierr);
32760f0b76eSHong Zhang       ierr = TSSetIJacobian(ts,appctx.A,appctx.A,IJacobianShell,&appctx);CHKERRQ(ierr);
32860f0b76eSHong Zhang     }
32960f0b76eSHong Zhang   }
33060f0b76eSHong Zhang 
33160f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33260f0b76eSHong Zhang      Set initial conditions
33360f0b76eSHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
33460f0b76eSHong Zhang   ierr = InitialConditions(da,x);CHKERRQ(ierr);
33560f0b76eSHong Zhang   ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
33660f0b76eSHong Zhang 
33760f0b76eSHong Zhang   /*
33860f0b76eSHong Zhang     Have the TS save its trajectory so that TSAdjointSolve() may be used
33960f0b76eSHong Zhang   */
34060f0b76eSHong Zhang   if (!forwardonly) { ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); }
34160f0b76eSHong Zhang 
34260f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34360f0b76eSHong Zhang      Set solver options
34460f0b76eSHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
34560f0b76eSHong Zhang   ierr = TSSetMaxTime(ts,200.0);CHKERRQ(ierr);
34660f0b76eSHong Zhang   ierr = TSSetTimeStep(ts,0.5);CHKERRQ(ierr);
34760f0b76eSHong Zhang   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
34860f0b76eSHong Zhang   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
34960f0b76eSHong Zhang 
35060f0b76eSHong Zhang   ierr = PetscTime(&v1);CHKERRQ(ierr);
35160f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35260f0b76eSHong Zhang      Solve ODE system
35360f0b76eSHong Zhang      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
35460f0b76eSHong Zhang   ierr = TSSolve(ts,x);CHKERRQ(ierr);
35560f0b76eSHong Zhang   if (!forwardonly) {
35660f0b76eSHong Zhang     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35760f0b76eSHong Zhang        Start the Adjoint model
35860f0b76eSHong Zhang        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
35960f0b76eSHong Zhang     ierr = VecDuplicate(x,&lambda[0]);CHKERRQ(ierr);
36060f0b76eSHong Zhang     /*   Reset initial conditions for the adjoint integration */
36160f0b76eSHong Zhang     ierr = InitializeLambda(da,lambda[0],0.5,0.5);CHKERRQ(ierr);
36260f0b76eSHong Zhang     ierr = TSSetCostGradients(ts,1,lambda,NULL);CHKERRQ(ierr);
36360f0b76eSHong Zhang     PetscLogStagePush(stage);
36460f0b76eSHong Zhang     ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
36560f0b76eSHong Zhang     PetscLogStagePop();
36660f0b76eSHong Zhang     ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr);
36760f0b76eSHong Zhang   }
36860f0b76eSHong Zhang   ierr = PetscTime(&v2);CHKERRQ(ierr);
36960f0b76eSHong Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD," %.3lf ",v2-v1);
37060f0b76eSHong Zhang 
37160f0b76eSHong Zhang   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37260f0b76eSHong Zhang      Free work space.  All PETSc objects should be destroyed when they
37360f0b76eSHong Zhang      are no longer needed.
37460f0b76eSHong Zhang    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
37560f0b76eSHong Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
37660f0b76eSHong Zhang   ierr = TSDestroy(&ts);CHKERRQ(ierr);
37760f0b76eSHong Zhang   ierr = DMDestroy(&da);CHKERRQ(ierr);
37860f0b76eSHong Zhang   if (mf) {
37960f0b76eSHong Zhang     ierr = VecDestroy(&mctx.U);CHKERRQ(ierr);
38060f0b76eSHong Zhang     /* ierr = VecDestroy(&mctx.Udot);CHKERRQ(ierr);*/
38160f0b76eSHong Zhang     ierr = MatDestroy(&appctx.A);CHKERRQ(ierr);
38260f0b76eSHong Zhang   }
38360f0b76eSHong Zhang   ierr = PetscFinalize();
38460f0b76eSHong Zhang   return ierr;
38560f0b76eSHong Zhang }
38660f0b76eSHong Zhang 
38760f0b76eSHong Zhang /* ------------------------------------------------------------------- */
38860f0b76eSHong Zhang PetscErrorCode RHSJacobianShell(TS ts,PetscReal t,Vec U,Mat A,Mat BB,void *ctx)
38960f0b76eSHong Zhang {
39060f0b76eSHong Zhang   MCtx           *mctx;
39160f0b76eSHong Zhang   PetscErrorCode ierr;
39260f0b76eSHong Zhang 
39360f0b76eSHong Zhang   PetscFunctionBegin;
39460f0b76eSHong Zhang   ierr = MatShellGetContext(A,&mctx);CHKERRQ(ierr);
39560f0b76eSHong Zhang   ierr = VecCopy(U,mctx->U);CHKERRQ(ierr);
39660f0b76eSHong Zhang   PetscFunctionReturn(0);
39760f0b76eSHong Zhang }
39860f0b76eSHong Zhang 
39960f0b76eSHong Zhang PetscErrorCode IJacobianShell(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat BB,void *ctx)
40060f0b76eSHong Zhang {
40160f0b76eSHong Zhang   MCtx           *mctx;
40260f0b76eSHong Zhang   PetscErrorCode ierr;
40360f0b76eSHong Zhang 
40460f0b76eSHong Zhang   PetscFunctionBegin;
40560f0b76eSHong Zhang   ierr = MatShellGetContext(A,&mctx);CHKERRQ(ierr);
40660f0b76eSHong Zhang   ierr = VecCopy(U,mctx->U);CHKERRQ(ierr);
40760f0b76eSHong Zhang   /* ierr = VecCopy(Udot,mctx->Udot);CHKERRQ(ierr); */
40860f0b76eSHong Zhang   mctx->shift = a;
40960f0b76eSHong Zhang   PetscFunctionReturn(0);
41060f0b76eSHong Zhang }
41160f0b76eSHong Zhang 
41260f0b76eSHong Zhang /* ------------------------------------------------------------------- */
41360f0b76eSHong Zhang PetscErrorCode InitialConditions(DM da,Vec U)
41460f0b76eSHong Zhang {
41560f0b76eSHong Zhang   PetscErrorCode ierr;
41660f0b76eSHong Zhang   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
41760f0b76eSHong Zhang   Field          **u;
41860f0b76eSHong Zhang   PetscReal      hx,hy,x,y;
41960f0b76eSHong Zhang 
42060f0b76eSHong Zhang   PetscFunctionBegin;
42160f0b76eSHong Zhang   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
42260f0b76eSHong Zhang 
42360f0b76eSHong Zhang   hx = 2.5/(PetscReal)Mx;
42460f0b76eSHong Zhang   hy = 2.5/(PetscReal)My;
42560f0b76eSHong Zhang 
42660f0b76eSHong Zhang   /*
42760f0b76eSHong Zhang      Get pointers to vector data
42860f0b76eSHong Zhang   */
42960f0b76eSHong Zhang   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
43060f0b76eSHong Zhang 
43160f0b76eSHong Zhang   /*
43260f0b76eSHong Zhang      Get local grid boundaries
43360f0b76eSHong Zhang   */
43460f0b76eSHong Zhang   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
43560f0b76eSHong Zhang 
43660f0b76eSHong Zhang   /*
43760f0b76eSHong Zhang      Compute function over the locally owned part of the grid
43860f0b76eSHong Zhang   */
43960f0b76eSHong Zhang   for (j=ys; j<ys+ym; j++) {
44060f0b76eSHong Zhang     y = j*hy;
44160f0b76eSHong Zhang     for (i=xs; i<xs+xm; i++) {
44260f0b76eSHong Zhang       x = i*hx;
44360f0b76eSHong 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);
44460f0b76eSHong Zhang       else u[j][i].v = 0.0;
44560f0b76eSHong Zhang 
44660f0b76eSHong Zhang       u[j][i].u = 1.0 - 2.0*u[j][i].v;
44760f0b76eSHong Zhang     }
44860f0b76eSHong Zhang   }
44960f0b76eSHong Zhang 
45060f0b76eSHong Zhang   /*
45160f0b76eSHong Zhang      Restore vectors
45260f0b76eSHong Zhang   */
45360f0b76eSHong Zhang   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
45460f0b76eSHong Zhang   PetscFunctionReturn(0);
45560f0b76eSHong Zhang }
45660f0b76eSHong Zhang 
45760f0b76eSHong Zhang /*TEST
45860f0b76eSHong Zhang 
45960f0b76eSHong Zhang    build:
46060f0b76eSHong Zhang       depends: reaction_diffusion.c
46160f0b76eSHong Zhang       requires: !complex !single
46260f0b76eSHong Zhang 
46360f0b76eSHong Zhang    test:
46460f0b76eSHong Zhang       requires: cams
465*533251d3SHong 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
46660f0b76eSHong Zhang TEST*/
467