xref: /petsc/src/ts/tutorials/autodiff/adr_ex5adj_mf.cxx (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Demonstrates automatic, matrix-free Jacobian generation using ADOL-C for a time-dependent PDE in 2d, solved using implicit timestepping.\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown /*
4*c4762a1bSJed Brown    Concepts: TS^time-dependent nonlinear problems
5*c4762a1bSJed Brown    Concepts: Automatic differentiation using ADOL-C
6*c4762a1bSJed Brown    Concepts: Matrix-free methods
7*c4762a1bSJed Brown */
8*c4762a1bSJed Brown /*
9*c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-adolc.
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown    For documentation on ADOL-C, see
12*c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
13*c4762a1bSJed Brown */
14*c4762a1bSJed Brown /* ------------------------------------------------------------------------
15*c4762a1bSJed Brown   See ../advection-diffusion-reaction/ex5 for a description of the problem
16*c4762a1bSJed Brown   ------------------------------------------------------------------------- */
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown #include <petscdmda.h>
19*c4762a1bSJed Brown #include <petscts.h>
20*c4762a1bSJed Brown #include "adolc-utils/init.cxx"
21*c4762a1bSJed Brown #include "adolc-utils/matfree.cxx"
22*c4762a1bSJed Brown #include <adolc/adolc.h>
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown /* (Passive) field for the two variables */
25*c4762a1bSJed Brown typedef struct {
26*c4762a1bSJed Brown   PetscScalar u,v;
27*c4762a1bSJed Brown } Field;
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown /* Active field for the two variables */
30*c4762a1bSJed Brown typedef struct {
31*c4762a1bSJed Brown   adouble u,v;
32*c4762a1bSJed Brown } AField;
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown /* Application context */
35*c4762a1bSJed Brown typedef struct {
36*c4762a1bSJed Brown   PetscReal D1,D2,gamma,kappa;
37*c4762a1bSJed Brown   AField    **u_a,**f_a;
38*c4762a1bSJed Brown   AdolcCtx  *adctx; /* Automatic differentation support */
39*c4762a1bSJed Brown } AppCtx;
40*c4762a1bSJed Brown 
41*c4762a1bSJed Brown extern PetscErrorCode InitialConditions(DM da,Vec U);
42*c4762a1bSJed Brown extern PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y);
43*c4762a1bSJed Brown extern PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info,PetscReal t,Field**u,Field**udot,Field**f,void *ptr);
44*c4762a1bSJed Brown extern PetscErrorCode IFunctionActive(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr);
45*c4762a1bSJed Brown extern PetscErrorCode IJacobianMatFree(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A_shell,Mat B,void *ctx);
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown int main(int argc,char **argv)
48*c4762a1bSJed Brown {
49*c4762a1bSJed Brown   TS             ts;                  /* ODE integrator */
50*c4762a1bSJed Brown   Vec            x,r;                 /* solution, residual */
51*c4762a1bSJed Brown   PetscErrorCode ierr;
52*c4762a1bSJed Brown   DM             da;
53*c4762a1bSJed Brown   AppCtx         appctx;              /* Application context */
54*c4762a1bSJed Brown   AdolcMatCtx    matctx;              /* Matrix (free) context */
55*c4762a1bSJed Brown   Vec            lambda[1];
56*c4762a1bSJed Brown   PetscBool      forwardonly=PETSC_FALSE;
57*c4762a1bSJed Brown   Mat            A;                   /* (Matrix free) Jacobian matrix */
58*c4762a1bSJed Brown   PetscInt       gxm,gym;
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61*c4762a1bSJed Brown      Initialize program
62*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
64*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL);CHKERRQ(ierr);
65*c4762a1bSJed Brown   PetscFunctionBeginUser;
66*c4762a1bSJed Brown   appctx.D1    = 8.0e-5;
67*c4762a1bSJed Brown   appctx.D2    = 4.0e-5;
68*c4762a1bSJed Brown   appctx.gamma = .024;
69*c4762a1bSJed Brown   appctx.kappa = .06;
70*c4762a1bSJed Brown   ierr = PetscLogEventRegister("df/dx forward",MAT_CLASSID,&matctx.event1);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = PetscLogEventRegister("df/d(xdot) forward",MAT_CLASSID,&matctx.event2);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = PetscLogEventRegister("df/dx reverse",MAT_CLASSID,&matctx.event3);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = PetscLogEventRegister("df/d(xdot) reverse",MAT_CLASSID,&matctx.event4);CHKERRQ(ierr);
74*c4762a1bSJed Brown 
75*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76*c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
77*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78*c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,65,65,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = DMDASetFieldName(da,0,"u");CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = DMDASetFieldName(da,1,"v");CHKERRQ(ierr);
83*c4762a1bSJed Brown 
84*c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85*c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
86*c4762a1bSJed Brown      vectors that are the same types
87*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
90*c4762a1bSJed Brown 
91*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92*c4762a1bSJed Brown     Create matrix free context and specify usage of PETSc-ADOL-C drivers
93*c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94*c4762a1bSJed Brown   ierr = DMSetMatType(da,MATSHELL);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = DMCreateMatrix(da,&A);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = MatShellSetContext(A,&matctx);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = MatShellSetOperation(A,MATOP_MULT,(void (*)(void))PetscAdolcIJacobianVectorProductIDMass);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void (*)(void))PetscAdolcIJacobianTransposeVectorProductIDMass);CHKERRQ(ierr);
99*c4762a1bSJed Brown   ierr = VecDuplicate(x,&matctx.X);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = VecDuplicate(x,&matctx.Xdot);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&matctx.localX0);CHKERRQ(ierr);
102*c4762a1bSJed Brown 
103*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104*c4762a1bSJed Brown      Create timestepping solver context
105*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = DMDATSSetIFunctionLocal(da,INSERT_VALUES,(DMDATSIFunctionLocal)IFunctionLocalPassive,&appctx);CHKERRQ(ierr);
111*c4762a1bSJed Brown 
112*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113*c4762a1bSJed Brown     Some data required for matrix-free context
114*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115*c4762a1bSJed Brown   ierr = DMDAGetGhostCorners(da,NULL,NULL,NULL,&gxm,&gym,NULL);CHKERRQ(ierr);
116*c4762a1bSJed Brown   matctx.m = 2*gxm*gym;matctx.n = 2*gxm*gym; /* Number of dependent and independent variables */
117*c4762a1bSJed Brown   matctx.flg = PETSC_FALSE;                  /* Flag for reverse mode */
118*c4762a1bSJed Brown   matctx.tag1 = 1;                           /* Tape identifier */
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121*c4762a1bSJed Brown      Trace function just once
122*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123*c4762a1bSJed Brown   ierr = PetscNew(&appctx.adctx);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = IFunctionActive(ts,1.,x,matctx.Xdot,r,&appctx);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = PetscFree(appctx.adctx);CHKERRQ(ierr);
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128*c4762a1bSJed Brown      Set Jacobian. In this case, IJacobian simply acts to pass context
129*c4762a1bSJed Brown      information to the matrix-free Jacobian vector product.
130*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131*c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,A,A,IJacobianMatFree,&appctx);CHKERRQ(ierr);
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134*c4762a1bSJed Brown      Set initial conditions
135*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136*c4762a1bSJed Brown   ierr = InitialConditions(da,x);CHKERRQ(ierr);
137*c4762a1bSJed Brown   ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140*c4762a1bSJed Brown     Have the TS save its trajectory so that TSAdjointSolve() may be used
141*c4762a1bSJed Brown     and set solver options
142*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143*c4762a1bSJed Brown   if (!forwardonly) {
144*c4762a1bSJed Brown     ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
145*c4762a1bSJed Brown     ierr = TSSetMaxTime(ts,200.0);CHKERRQ(ierr);
146*c4762a1bSJed Brown     ierr = TSSetTimeStep(ts,0.5);CHKERRQ(ierr);
147*c4762a1bSJed Brown   } else {
148*c4762a1bSJed Brown     ierr = TSSetMaxTime(ts,2000.0);CHKERRQ(ierr);
149*c4762a1bSJed Brown     ierr = TSSetTimeStep(ts,10);CHKERRQ(ierr);
150*c4762a1bSJed Brown   }
151*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155*c4762a1bSJed Brown      Solve ODE system
156*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157*c4762a1bSJed Brown   ierr = TSSolve(ts,x);CHKERRQ(ierr);
158*c4762a1bSJed Brown   if (!forwardonly) {
159*c4762a1bSJed Brown     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160*c4762a1bSJed Brown        Start the Adjoint model
161*c4762a1bSJed Brown        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
162*c4762a1bSJed Brown     ierr = VecDuplicate(x,&lambda[0]);CHKERRQ(ierr);
163*c4762a1bSJed Brown     /*   Reset initial conditions for the adjoint integration */
164*c4762a1bSJed Brown     ierr = InitializeLambda(da,lambda[0],0.5,0.5);CHKERRQ(ierr);
165*c4762a1bSJed Brown     ierr = TSSetCostGradients(ts,1,lambda,NULL);CHKERRQ(ierr);
166*c4762a1bSJed Brown     ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
167*c4762a1bSJed Brown     ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr);
168*c4762a1bSJed Brown   }
169*c4762a1bSJed Brown 
170*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
172*c4762a1bSJed Brown      are no longer needed.
173*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
174*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&matctx.localX0);CHKERRQ(ierr);
175*c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
176*c4762a1bSJed Brown   ierr = VecDestroy(&matctx.X);CHKERRQ(ierr);
177*c4762a1bSJed Brown   ierr = VecDestroy(&matctx.Xdot);CHKERRQ(ierr);
178*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
179*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
180*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
181*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
182*c4762a1bSJed Brown 
183*c4762a1bSJed Brown   ierr = PetscFinalize();
184*c4762a1bSJed Brown   return ierr;
185*c4762a1bSJed Brown }
186*c4762a1bSJed Brown 
187*c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U)
188*c4762a1bSJed Brown {
189*c4762a1bSJed Brown   PetscErrorCode ierr;
190*c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
191*c4762a1bSJed Brown   Field          **u;
192*c4762a1bSJed Brown   PetscReal      hx,hy,x,y;
193*c4762a1bSJed Brown 
194*c4762a1bSJed Brown   PetscFunctionBegin;
195*c4762a1bSJed Brown   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);
196*c4762a1bSJed Brown 
197*c4762a1bSJed Brown   hx = 2.5/(PetscReal)Mx;
198*c4762a1bSJed Brown   hy = 2.5/(PetscReal)My;
199*c4762a1bSJed Brown 
200*c4762a1bSJed Brown   /*
201*c4762a1bSJed Brown      Get pointers to vector data
202*c4762a1bSJed Brown   */
203*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,U,&u);CHKERRQ(ierr);
204*c4762a1bSJed Brown 
205*c4762a1bSJed Brown   /*
206*c4762a1bSJed Brown      Get local grid boundaries
207*c4762a1bSJed Brown   */
208*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
209*c4762a1bSJed Brown 
210*c4762a1bSJed Brown   /*
211*c4762a1bSJed Brown      Compute function over the locally owned part of the grid
212*c4762a1bSJed Brown   */
213*c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
214*c4762a1bSJed Brown     y = j*hy;
215*c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
216*c4762a1bSJed Brown       x = i*hx;
217*c4762a1bSJed Brown       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);
218*c4762a1bSJed Brown       else u[j][i].v = 0.0;
219*c4762a1bSJed Brown 
220*c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0*u[j][i].v;
221*c4762a1bSJed Brown     }
222*c4762a1bSJed Brown   }
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown   /*
225*c4762a1bSJed Brown      Restore vectors
226*c4762a1bSJed Brown   */
227*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,U,&u);CHKERRQ(ierr);
228*c4762a1bSJed Brown   PetscFunctionReturn(0);
229*c4762a1bSJed Brown }
230*c4762a1bSJed Brown 
231*c4762a1bSJed Brown PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y)
232*c4762a1bSJed Brown {
233*c4762a1bSJed Brown    PetscInt i,j,Mx,My,xs,ys,xm,ym;
234*c4762a1bSJed Brown    PetscErrorCode ierr;
235*c4762a1bSJed Brown    Field **l;
236*c4762a1bSJed Brown    PetscFunctionBegin;
237*c4762a1bSJed Brown 
238*c4762a1bSJed Brown    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);
239*c4762a1bSJed Brown    /* locate the global i index for x and j index for y */
240*c4762a1bSJed Brown    i = (PetscInt)(x*(Mx-1));
241*c4762a1bSJed Brown    j = (PetscInt)(y*(My-1));
242*c4762a1bSJed Brown    ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
243*c4762a1bSJed Brown 
244*c4762a1bSJed Brown    if (xs <= i && i < xs+xm && ys <= j && j < ys+ym) {
245*c4762a1bSJed Brown      /* the i,j vertex is on this process */
246*c4762a1bSJed Brown      ierr = DMDAVecGetArray(da,lambda,&l);CHKERRQ(ierr);
247*c4762a1bSJed Brown      l[j][i].u = 1.0;
248*c4762a1bSJed Brown      l[j][i].v = 1.0;
249*c4762a1bSJed Brown      ierr = DMDAVecRestoreArray(da,lambda,&l);CHKERRQ(ierr);
250*c4762a1bSJed Brown    }
251*c4762a1bSJed Brown    PetscFunctionReturn(0);
252*c4762a1bSJed Brown }
253*c4762a1bSJed Brown 
254*c4762a1bSJed Brown PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info,PetscReal t,Field**u,Field**udot,Field**f,void *ptr)
255*c4762a1bSJed Brown {
256*c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
257*c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym;
258*c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
259*c4762a1bSJed Brown   PetscScalar    uc,uxx,uyy,vc,vxx,vyy;
260*c4762a1bSJed Brown   PetscErrorCode ierr;
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown   PetscFunctionBegin;
263*c4762a1bSJed Brown   hx = 2.50/(PetscReal)(info->mx); sx = 1.0/(hx*hx);
264*c4762a1bSJed Brown   hy = 2.50/(PetscReal)(info->my); sy = 1.0/(hy*hy);
265*c4762a1bSJed Brown 
266*c4762a1bSJed Brown   /* Get local grid boundaries */
267*c4762a1bSJed Brown   xs = info->xs; xm = info->xm; ys = info->ys; ym = info->ym;
268*c4762a1bSJed Brown 
269*c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
270*c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
271*c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
272*c4762a1bSJed Brown       uc        = u[j][i].u;
273*c4762a1bSJed Brown       uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
274*c4762a1bSJed Brown       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
275*c4762a1bSJed Brown       vc        = u[j][i].v;
276*c4762a1bSJed Brown       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
277*c4762a1bSJed Brown       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
278*c4762a1bSJed Brown       f[j][i].u = udot[j][i].u - appctx->D1*(uxx + uyy) + uc*vc*vc - appctx->gamma*(1.0 - uc);
279*c4762a1bSJed Brown       f[j][i].v = udot[j][i].v - appctx->D2*(vxx + vyy) - uc*vc*vc + (appctx->gamma + appctx->kappa)*vc;
280*c4762a1bSJed Brown     }
281*c4762a1bSJed Brown   }
282*c4762a1bSJed Brown   ierr = PetscLogFlops(16*xm*ym);CHKERRQ(ierr);
283*c4762a1bSJed Brown   PetscFunctionReturn(0);
284*c4762a1bSJed Brown }
285*c4762a1bSJed Brown 
286*c4762a1bSJed Brown PetscErrorCode IFunctionActive(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
287*c4762a1bSJed Brown {
288*c4762a1bSJed Brown   PetscErrorCode ierr;
289*c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
290*c4762a1bSJed Brown   DM             da;
291*c4762a1bSJed Brown   DMDALocalInfo  info;
292*c4762a1bSJed Brown   Field          **u,**f,**udot;
293*c4762a1bSJed Brown   Vec            localU;
294*c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,gxs,gys,gxm,gym;
295*c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
296*c4762a1bSJed Brown   adouble        uc,uxx,uyy,vc,vxx,vyy;
297*c4762a1bSJed Brown   AField         **f_a,*f_c,**u_a,*u_c;
298*c4762a1bSJed Brown   PetscScalar    dummy;
299*c4762a1bSJed Brown 
300*c4762a1bSJed Brown   PetscFunctionBegin;
301*c4762a1bSJed Brown 
302*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
303*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
304*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
305*c4762a1bSJed Brown   hx = 2.50/(PetscReal)(info.mx); sx = 1.0/(hx*hx);
306*c4762a1bSJed Brown   hy = 2.50/(PetscReal)(info.my); sy = 1.0/(hy*hy);
307*c4762a1bSJed Brown   xs = info.xs; xm = info.xm; gxs = info.gxs; gxm = info.gxm;
308*c4762a1bSJed Brown   ys = info.ys; ym = info.ym; gys = info.gys; gym = info.gym;
309*c4762a1bSJed Brown 
310*c4762a1bSJed Brown   /*
311*c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
312*c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
313*c4762a1bSJed Brown      By placing code between these two statements, computations can be
314*c4762a1bSJed Brown      done while messages are in transition.
315*c4762a1bSJed Brown   */
316*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
317*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
318*c4762a1bSJed Brown 
319*c4762a1bSJed Brown   /*
320*c4762a1bSJed Brown      Get pointers to vector data
321*c4762a1bSJed Brown   */
322*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
323*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,F,&f);CHKERRQ(ierr);
324*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Udot,&udot);CHKERRQ(ierr);
325*c4762a1bSJed Brown 
326*c4762a1bSJed Brown   /*
327*c4762a1bSJed Brown     Create contiguous 1-arrays of AFields
328*c4762a1bSJed Brown 
329*c4762a1bSJed Brown     NOTE: Memory for ADOL-C active variables (such as adouble and AField)
330*c4762a1bSJed Brown           cannot be allocated using PetscMalloc, as this does not call the
331*c4762a1bSJed Brown           relevant class constructor. Instead, we use the C++ keyword `new`.
332*c4762a1bSJed Brown   */
333*c4762a1bSJed Brown   u_c = new AField[info.gxm*info.gym];
334*c4762a1bSJed Brown   f_c = new AField[info.gxm*info.gym];
335*c4762a1bSJed Brown 
336*c4762a1bSJed Brown   /* Create corresponding 2-arrays of AFields */
337*c4762a1bSJed Brown   u_a = new AField*[info.gym];
338*c4762a1bSJed Brown   f_a = new AField*[info.gym];
339*c4762a1bSJed Brown 
340*c4762a1bSJed Brown   /* Align indices between array types to endow 2d array with ghost points */
341*c4762a1bSJed Brown   ierr = GiveGhostPoints(da,u_c,&u_a);CHKERRQ(ierr);
342*c4762a1bSJed Brown   ierr = GiveGhostPoints(da,f_c,&f_a);CHKERRQ(ierr);
343*c4762a1bSJed Brown 
344*c4762a1bSJed Brown   trace_on(1);  /* Start of active section on tape 1 */
345*c4762a1bSJed Brown 
346*c4762a1bSJed Brown   /*
347*c4762a1bSJed Brown     Mark independence
348*c4762a1bSJed Brown 
349*c4762a1bSJed Brown     NOTE: Ghost points are marked as independent, in place of the points they represent on
350*c4762a1bSJed Brown           other processors / on other boundaries.
351*c4762a1bSJed Brown   */
352*c4762a1bSJed Brown   for (j=gys; j<gys+gym; j++) {
353*c4762a1bSJed Brown     for (i=gxs; i<gxs+gxm; i++) {
354*c4762a1bSJed Brown       u_a[j][i].u <<= u[j][i].u;
355*c4762a1bSJed Brown       u_a[j][i].v <<= u[j][i].v;
356*c4762a1bSJed Brown     }
357*c4762a1bSJed Brown   }
358*c4762a1bSJed Brown 
359*c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
360*c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
361*c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
362*c4762a1bSJed Brown       uc        = u_a[j][i].u;
363*c4762a1bSJed Brown       uxx       = (-2.0*uc + u_a[j][i-1].u + u_a[j][i+1].u)*sx;
364*c4762a1bSJed Brown       uyy       = (-2.0*uc + u_a[j-1][i].u + u_a[j+1][i].u)*sy;
365*c4762a1bSJed Brown       vc        = u_a[j][i].v;
366*c4762a1bSJed Brown       vxx       = (-2.0*vc + u_a[j][i-1].v + u_a[j][i+1].v)*sx;
367*c4762a1bSJed Brown       vyy       = (-2.0*vc + u_a[j-1][i].v + u_a[j+1][i].v)*sy;
368*c4762a1bSJed Brown       f_a[j][i].u = udot[j][i].u - appctx->D1*(uxx + uyy) + uc*vc*vc - appctx->gamma*(1.0 - uc);
369*c4762a1bSJed Brown       f_a[j][i].v = udot[j][i].v - appctx->D2*(vxx + vyy) - uc*vc*vc + (appctx->gamma + appctx->kappa)*vc;
370*c4762a1bSJed Brown     }
371*c4762a1bSJed Brown   }
372*c4762a1bSJed Brown 
373*c4762a1bSJed Brown   /*
374*c4762a1bSJed Brown     Mark dependence
375*c4762a1bSJed Brown 
376*c4762a1bSJed Brown     NOTE: Marking dependence of dummy variables makes the index notation much simpler when forming
377*c4762a1bSJed Brown           the Jacobian later.
378*c4762a1bSJed Brown   */
379*c4762a1bSJed Brown   for (j=gys; j<gys+gym; j++) {
380*c4762a1bSJed Brown     for (i=gxs; i<gxs+gxm; i++) {
381*c4762a1bSJed Brown       if ((i < xs) || (i >= xs+xm) || (j < ys) || (j >= ys+ym)) {
382*c4762a1bSJed Brown         f_a[j][i].u >>= dummy;
383*c4762a1bSJed Brown         f_a[j][i].v >>= dummy;
384*c4762a1bSJed Brown       } else {
385*c4762a1bSJed Brown         f_a[j][i].u >>= f[j][i].u;
386*c4762a1bSJed Brown         f_a[j][i].v >>= f[j][i].v;
387*c4762a1bSJed Brown       }
388*c4762a1bSJed Brown     }
389*c4762a1bSJed Brown   }
390*c4762a1bSJed Brown   trace_off();  /* End of active section */
391*c4762a1bSJed Brown   ierr = PetscLogFlops(16*xm*ym);CHKERRQ(ierr);
392*c4762a1bSJed Brown 
393*c4762a1bSJed Brown   /* Restore vectors */
394*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,F,&f);CHKERRQ(ierr);
395*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
396*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Udot,&udot);CHKERRQ(ierr);
397*c4762a1bSJed Brown 
398*c4762a1bSJed Brown   /* Destroy AFields appropriately */
399*c4762a1bSJed Brown   f_a += info.gys;
400*c4762a1bSJed Brown   u_a += info.gys;
401*c4762a1bSJed Brown   delete[] f_a;
402*c4762a1bSJed Brown   delete[] u_a;
403*c4762a1bSJed Brown   delete[] f_c;
404*c4762a1bSJed Brown   delete[] u_c;
405*c4762a1bSJed Brown 
406*c4762a1bSJed Brown   PetscFunctionReturn(0);
407*c4762a1bSJed Brown }
408*c4762a1bSJed Brown 
409*c4762a1bSJed Brown /*
410*c4762a1bSJed Brown   Simply acts to pass TS information to the AdolcMatCtx
411*c4762a1bSJed Brown */
412*c4762a1bSJed Brown PetscErrorCode IJacobianMatFree(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A_shell,Mat B,void *ctx)
413*c4762a1bSJed Brown {
414*c4762a1bSJed Brown   AdolcMatCtx       *mctx;
415*c4762a1bSJed Brown   PetscErrorCode    ierr;
416*c4762a1bSJed Brown   DM                da;
417*c4762a1bSJed Brown 
418*c4762a1bSJed Brown   PetscFunctionBeginUser;
419*c4762a1bSJed Brown   ierr = MatShellGetContext(A_shell,(void **)&mctx);CHKERRQ(ierr);
420*c4762a1bSJed Brown 
421*c4762a1bSJed Brown   mctx->time  = t;
422*c4762a1bSJed Brown   mctx->shift = a;
423*c4762a1bSJed Brown   if (mctx->ts != ts) mctx->ts = ts;
424*c4762a1bSJed Brown   ierr = VecCopy(X,mctx->X);CHKERRQ(ierr);
425*c4762a1bSJed Brown   ierr = VecCopy(Xdot,mctx->Xdot);CHKERRQ(ierr);
426*c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
427*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,mctx->X,INSERT_VALUES,mctx->localX0);CHKERRQ(ierr);
428*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,mctx->X,INSERT_VALUES,mctx->localX0);CHKERRQ(ierr);
429*c4762a1bSJed Brown   PetscFunctionReturn(0);
430*c4762a1bSJed Brown }
431*c4762a1bSJed Brown 
432*c4762a1bSJed Brown /*TEST
433*c4762a1bSJed Brown 
434*c4762a1bSJed Brown   build:
435*c4762a1bSJed Brown     requires: double !complex adolc
436*c4762a1bSJed Brown 
437*c4762a1bSJed Brown   test:
438*c4762a1bSJed Brown     suffix: 1
439*c4762a1bSJed Brown     args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian
440*c4762a1bSJed Brown     output_file: output/adr_ex5adj_mf_1.out
441*c4762a1bSJed Brown 
442*c4762a1bSJed Brown   test:
443*c4762a1bSJed Brown     suffix: 2
444*c4762a1bSJed Brown     nsize: 4
445*c4762a1bSJed Brown     args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor
446*c4762a1bSJed Brown     output_file: output/adr_ex5adj_mf_2.out
447*c4762a1bSJed Brown 
448*c4762a1bSJed Brown TEST*/
449