xref: /petsc/src/ts/tutorials/autodiff/adr_ex5adj.cxx (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static char help[] = "Demonstrates adjoint sensitivity analysis for Reaction-Diffusion Equations.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-adolc.
5c4762a1bSJed Brown 
6c4762a1bSJed Brown    For documentation on ADOL-C, see
7c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
8c4762a1bSJed Brown 
9c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-colpack
10c4762a1bSJed Brown 
11c4762a1bSJed Brown    For documentation on ColPack, see
12c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/git.colpack/README.md
13c4762a1bSJed Brown */
14c4762a1bSJed Brown /* ------------------------------------------------------------------------
15c4762a1bSJed Brown   See ../advection-diffusion-reaction/ex5 for a description of the problem
16c4762a1bSJed Brown   ------------------------------------------------------------------------- */
17c4762a1bSJed Brown 
18c4762a1bSJed Brown /*
19c4762a1bSJed Brown   Runtime options:
20c4762a1bSJed Brown 
21c4762a1bSJed Brown     Solver:
22c4762a1bSJed Brown       -forwardonly       - Run the forward simulation without adjoint.
23c4762a1bSJed Brown       -implicitform      - Provide IFunction and IJacobian to TS, if not set, RHSFunction and RHSJacobian will be used.
24c4762a1bSJed Brown       -aijpc             - Set the preconditioner matrix to be aij (the Jacobian matrix can be of a different type such as ELL).
25c4762a1bSJed Brown 
26c4762a1bSJed Brown     Jacobian generation:
27c4762a1bSJed Brown       -jacobian_by_hand  - Use the hand-coded Jacobian of ex5.c, rather than generating it automatically.
28c4762a1bSJed Brown       -no_annotation     - Do not annotate ADOL-C active variables. (Should be used alongside -jacobian_by_hand.)
29c4762a1bSJed Brown       -adolc_sparse      - Calculate Jacobian in compressed format, using ADOL-C sparse functionality.
30c4762a1bSJed Brown       -adolc_sparse_view - Print sparsity pattern used by -adolc_sparse option.
31c4762a1bSJed Brown  */
32c4762a1bSJed Brown /*
33c4762a1bSJed Brown   NOTE: If -adolc_sparse option is used, at least four processors should be used, so that no processor owns boundaries which are
34c4762a1bSJed Brown         identified by the periodic boundary conditions. The number of grid points in both x- and y-directions should be multiples
35c4762a1bSJed Brown         of 5, in order for the 5-point stencil to be cleanly parallelised.
36c4762a1bSJed Brown */
37c4762a1bSJed Brown 
38c4762a1bSJed Brown #include <petscdmda.h>
39c4762a1bSJed Brown #include <petscts.h>
40c4762a1bSJed Brown #include "adolc-utils/drivers.cxx"
41c4762a1bSJed Brown #include <adolc/adolc.h>
42c4762a1bSJed Brown 
43c4762a1bSJed Brown /* (Passive) field for the two variables */
44c4762a1bSJed Brown typedef struct {
45c4762a1bSJed Brown   PetscScalar u,v;
46c4762a1bSJed Brown } Field;
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /* Active field for the two variables */
49c4762a1bSJed Brown typedef struct {
50c4762a1bSJed Brown   adouble u,v;
51c4762a1bSJed Brown } AField;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown /* Application context */
54c4762a1bSJed Brown typedef struct {
55c4762a1bSJed Brown   PetscReal D1,D2,gamma,kappa;
56c4762a1bSJed Brown   AField    **u_a,**f_a;
57c4762a1bSJed Brown   PetscBool aijpc;
58c4762a1bSJed Brown   AdolcCtx  *adctx; /* Automatic differentation support */
59c4762a1bSJed Brown } AppCtx;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown extern PetscErrorCode InitialConditions(DM da,Vec U);
62c4762a1bSJed Brown extern PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y);
63c4762a1bSJed Brown extern PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info,PetscReal t,Field**u,Field**udot,Field**f,void *ptr);
64c4762a1bSJed Brown extern PetscErrorCode IFunctionActive(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr);
65c4762a1bSJed Brown extern PetscErrorCode RHSFunctionActive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr);
66c4762a1bSJed Brown extern PetscErrorCode RHSFunctionPassive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr);
67c4762a1bSJed Brown extern PetscErrorCode IJacobianByHand(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx);
68c4762a1bSJed Brown extern PetscErrorCode IJacobianAdolc(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx);
69c4762a1bSJed Brown extern PetscErrorCode RHSJacobianByHand(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx);
70c4762a1bSJed Brown extern PetscErrorCode RHSJacobianAdolc(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx);
71c4762a1bSJed Brown 
72c4762a1bSJed Brown int main(int argc,char **argv)
73c4762a1bSJed Brown {
74410585f6SHong Zhang   TS             ts;
75410585f6SHong Zhang   Vec            x,r,xdot;
76c4762a1bSJed Brown   DM             da;
77c4762a1bSJed Brown   AppCtx         appctx;
78c4762a1bSJed Brown   AdolcCtx       *adctx;
79c4762a1bSJed Brown   Vec            lambda[1];
80c4762a1bSJed Brown   PetscBool      forwardonly=PETSC_FALSE,implicitform=PETSC_FALSE,byhand=PETSC_FALSE;
81c4762a1bSJed Brown   PetscInt       gxm,gym,i,dofs = 2,ctrl[3] = {0,0,0};
82c4762a1bSJed Brown   PetscScalar    **Seed = NULL,**Rec = NULL,*u_vec;
83c4762a1bSJed Brown   unsigned int   **JP = NULL;
84c4762a1bSJed Brown   ISColoring     iscoloring;
85c4762a1bSJed Brown 
86*327415f7SBarry Smith   PetscFunctionBeginUser;
879566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
889566063dSJacob Faibussowitsch   PetscCall(PetscNew(&adctx));
899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL));
909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL));
91c4762a1bSJed Brown   appctx.aijpc = PETSC_FALSE,adctx->no_an = PETSC_FALSE,adctx->sparse = PETSC_FALSE,adctx->sparse_view = PETSC_FALSE;adctx->sparse_view_done = PETSC_FALSE;
929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-aijpc",&appctx.aijpc,NULL));
939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-no_annotation",&adctx->no_an,NULL));
949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-jacobian_by_hand",&byhand,NULL));
959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-adolc_sparse",&adctx->sparse,NULL));
969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-adolc_sparse_view",&adctx->sparse_view,NULL));
97c4762a1bSJed Brown   appctx.D1    = 8.0e-5;
98c4762a1bSJed Brown   appctx.D2    = 4.0e-5;
99c4762a1bSJed Brown   appctx.gamma = .024;
100c4762a1bSJed Brown   appctx.kappa = .06;
101c4762a1bSJed Brown   appctx.adctx = adctx;
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
105c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1069566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,65,65,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da));
1079566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1089566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1099566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da,0,"u"));
1109566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da,1,"v"));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
114c4762a1bSJed Brown      vectors that are the same types
115c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1169566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&x));
1179566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&r));
1189566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&xdot));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121c4762a1bSJed Brown      Create timestepping solver context
122c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1239566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
1249566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSCN));
1259566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts,da));
1269566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
127c4762a1bSJed Brown   if (!implicitform) {
1289566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&appctx));
129c4762a1bSJed Brown   } else {
1309566063dSJacob Faibussowitsch     PetscCall(DMDATSSetIFunctionLocal(da,INSERT_VALUES,(DMDATSIFunctionLocal)IFunctionLocalPassive,&appctx));
131c4762a1bSJed Brown   }
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   if (!adctx->no_an) {
1349566063dSJacob Faibussowitsch     PetscCall(DMDAGetGhostCorners(da,NULL,NULL,NULL,&gxm,&gym,NULL));
135c4762a1bSJed Brown     adctx->m = dofs*gxm*gym;adctx->n = dofs*gxm*gym; /* Number of dependent and independent variables */
136c4762a1bSJed Brown 
137c4762a1bSJed Brown     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138c4762a1bSJed Brown        Trace function(s) just once
139c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140c4762a1bSJed Brown     if (!implicitform) {
1419566063dSJacob Faibussowitsch       PetscCall(RHSFunctionActive(ts,1.0,x,r,&appctx));
142c4762a1bSJed Brown     } else {
1439566063dSJacob Faibussowitsch       PetscCall(IFunctionActive(ts,1.0,x,xdot,r,&appctx));
144c4762a1bSJed Brown     }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147c4762a1bSJed Brown       In the case where ADOL-C generates the Jacobian in compressed format,
148c4762a1bSJed Brown       seed and recovery matrices are required. Since the sparsity structure
149c4762a1bSJed Brown       of the Jacobian does not change over the course of the time
150c4762a1bSJed Brown       integration, we can save computational effort by only generating
151c4762a1bSJed Brown       these objects once.
152c4762a1bSJed Brown        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153c4762a1bSJed Brown     if (adctx->sparse) {
154c4762a1bSJed Brown       /*
155c4762a1bSJed Brown          Generate sparsity pattern
156c4762a1bSJed Brown 
157c4762a1bSJed Brown          One way to do this is to use the Jacobian pattern driver `jac_pat`
158c4762a1bSJed Brown          provided by ColPack. Otherwise, it is the responsibility of the user
159c4762a1bSJed Brown          to obtain the sparsity pattern.
160c4762a1bSJed Brown       */
1619566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(adctx->n,&u_vec));
162c4762a1bSJed Brown       JP = (unsigned int **) malloc(adctx->m*sizeof(unsigned int*));
163c4762a1bSJed Brown       jac_pat(1,adctx->m,adctx->n,u_vec,JP,ctrl);
164c4762a1bSJed Brown       if (adctx->sparse_view) {
1659566063dSJacob Faibussowitsch         PetscCall(PrintSparsity(MPI_COMM_WORLD,adctx->m,JP));
166c4762a1bSJed Brown       }
167c4762a1bSJed Brown 
168c4762a1bSJed Brown       /*
169c4762a1bSJed Brown         Extract a column colouring
170c4762a1bSJed Brown 
171c4762a1bSJed Brown         For problems using DMDA, colourings can be extracted directly, as
172c4762a1bSJed Brown         follows. There are also ColPack drivers available. Otherwise, it is the
173c4762a1bSJed Brown         responsibility of the user to obtain a suitable colouring.
174c4762a1bSJed Brown       */
1759566063dSJacob Faibussowitsch       PetscCall(DMCreateColoring(da,IS_COLORING_LOCAL,&iscoloring));
1769566063dSJacob Faibussowitsch       PetscCall(ISColoringGetIS(iscoloring,PETSC_USE_POINTER,&adctx->p,NULL));
177c4762a1bSJed Brown 
178c4762a1bSJed Brown       /* Generate seed matrix to propagate through the forward mode of AD */
1799566063dSJacob Faibussowitsch       PetscCall(AdolcMalloc2(adctx->n,adctx->p,&Seed));
1809566063dSJacob Faibussowitsch       PetscCall(GenerateSeedMatrix(iscoloring,Seed));
1819566063dSJacob Faibussowitsch       PetscCall(ISColoringDestroy(&iscoloring));
182c4762a1bSJed Brown 
183c4762a1bSJed Brown       /*
184c4762a1bSJed Brown         Generate recovery matrix, which is used to recover the Jacobian from
185c4762a1bSJed Brown         compressed format */
1869566063dSJacob Faibussowitsch       PetscCall(AdolcMalloc2(adctx->m,adctx->p,&Rec));
1879566063dSJacob Faibussowitsch       PetscCall(GetRecoveryMatrix(Seed,JP,adctx->m,adctx->p,Rec));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown       /* Store results and free workspace */
190c4762a1bSJed Brown       adctx->Rec = Rec;
191c4762a1bSJed Brown       for (i=0;i<adctx->m;i++)
192c4762a1bSJed Brown         free(JP[i]);
193c4762a1bSJed Brown       free(JP);
1949566063dSJacob Faibussowitsch       PetscCall(PetscFree(u_vec));
195c4762a1bSJed Brown 
196c4762a1bSJed Brown     } else {
197c4762a1bSJed Brown 
198c4762a1bSJed Brown       /*
199c4762a1bSJed Brown         In 'full' Jacobian mode, propagate an identity matrix through the
200c4762a1bSJed Brown         forward mode of AD.
201c4762a1bSJed Brown       */
202c4762a1bSJed Brown       adctx->p = adctx->n;
2039566063dSJacob Faibussowitsch       PetscCall(AdolcMalloc2(adctx->n,adctx->p,&Seed));
2049566063dSJacob Faibussowitsch       PetscCall(Identity(adctx->n,Seed));
205c4762a1bSJed Brown     }
206c4762a1bSJed Brown     adctx->Seed = Seed;
207c4762a1bSJed Brown   }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210c4762a1bSJed Brown      Set Jacobian
211c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212c4762a1bSJed Brown   if (!implicitform) {
213c4762a1bSJed Brown     if (!byhand) {
2149566063dSJacob Faibussowitsch       PetscCall(TSSetRHSJacobian(ts,NULL,NULL,RHSJacobianAdolc,&appctx));
215c4762a1bSJed Brown     } else {
2169566063dSJacob Faibussowitsch       PetscCall(TSSetRHSJacobian(ts,NULL,NULL,RHSJacobianByHand,&appctx));
217c4762a1bSJed Brown     }
218c4762a1bSJed Brown   } else {
219c4762a1bSJed Brown     if (appctx.aijpc) {
220c4762a1bSJed Brown       Mat A,B;
221c4762a1bSJed Brown 
2229566063dSJacob Faibussowitsch       PetscCall(DMSetMatType(da,MATSELL));
2239566063dSJacob Faibussowitsch       PetscCall(DMCreateMatrix(da,&A));
2249566063dSJacob Faibussowitsch       PetscCall(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B));
225c4762a1bSJed Brown       /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */
226c4762a1bSJed Brown       if (!byhand) {
2279566063dSJacob Faibussowitsch         PetscCall(TSSetIJacobian(ts,A,B,IJacobianAdolc,&appctx));
228c4762a1bSJed Brown       } else {
2299566063dSJacob Faibussowitsch         PetscCall(TSSetIJacobian(ts,A,B,IJacobianByHand,&appctx));
230c4762a1bSJed Brown       }
2319566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&A));
2329566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B));
233c4762a1bSJed Brown     } else {
234c4762a1bSJed Brown       if (!byhand) {
2359566063dSJacob Faibussowitsch         PetscCall(TSSetIJacobian(ts,NULL,NULL,IJacobianAdolc,&appctx));
236c4762a1bSJed Brown       } else {
2379566063dSJacob Faibussowitsch         PetscCall(TSSetIJacobian(ts,NULL,NULL,IJacobianByHand,&appctx));
238c4762a1bSJed Brown       }
239c4762a1bSJed Brown     }
240c4762a1bSJed Brown   }
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
243c4762a1bSJed Brown      Set initial conditions
244c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2459566063dSJacob Faibussowitsch   PetscCall(InitialConditions(da,x));
2469566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,x));
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
249c4762a1bSJed Brown     Have the TS save its trajectory so that TSAdjointSolve() may be used
250c4762a1bSJed Brown     and set solver options
251c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
252c4762a1bSJed Brown   if (!forwardonly) {
2539566063dSJacob Faibussowitsch     PetscCall(TSSetSaveTrajectory(ts));
2549566063dSJacob Faibussowitsch     PetscCall(TSSetMaxTime(ts,200.0));
2559566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts,0.5));
256c4762a1bSJed Brown   } else {
2579566063dSJacob Faibussowitsch     PetscCall(TSSetMaxTime(ts,2000.0));
2589566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts,10));
259c4762a1bSJed Brown   }
2609566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
2619566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
264c4762a1bSJed Brown      Solve ODE system
265c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2669566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,x));
267c4762a1bSJed Brown   if (!forwardonly) {
268c4762a1bSJed Brown     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269c4762a1bSJed Brown        Start the Adjoint model
270c4762a1bSJed Brown        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2719566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(x,&lambda[0]));
272c4762a1bSJed Brown     /*   Reset initial conditions for the adjoint integration */
2739566063dSJacob Faibussowitsch     PetscCall(InitializeLambda(da,lambda[0],0.5,0.5));
2749566063dSJacob Faibussowitsch     PetscCall(TSSetCostGradients(ts,1,lambda,NULL));
2759566063dSJacob Faibussowitsch     PetscCall(TSAdjointSolve(ts));
2769566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lambda[0]));
277c4762a1bSJed Brown   }
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
280c4762a1bSJed Brown      Free work space.
281c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&xdot));
2839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
2849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2859566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
286c4762a1bSJed Brown   if (!adctx->no_an) {
2879566063dSJacob Faibussowitsch     if (adctx->sparse) PetscCall(AdolcFree2(Rec));
2889566063dSJacob Faibussowitsch     PetscCall(AdolcFree2(Seed));
289c4762a1bSJed Brown   }
2909566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
2919566063dSJacob Faibussowitsch   PetscCall(PetscFree(adctx));
2929566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
293b122ec5aSJacob Faibussowitsch   return 0;
294c4762a1bSJed Brown }
295c4762a1bSJed Brown 
296c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U)
297c4762a1bSJed Brown {
298c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
299c4762a1bSJed Brown   Field          **u;
300c4762a1bSJed Brown   PetscReal      hx,hy,x,y;
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   PetscFunctionBegin;
3039566063dSJacob 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));
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   hx = 2.5/(PetscReal)Mx;
306c4762a1bSJed Brown   hy = 2.5/(PetscReal)My;
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   /*
309c4762a1bSJed Brown      Get pointers to vector data
310c4762a1bSJed Brown   */
3119566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,U,&u));
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   /*
314c4762a1bSJed Brown      Get local grid boundaries
315c4762a1bSJed Brown   */
3169566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   /*
319c4762a1bSJed Brown      Compute function over the locally owned part of the grid
320c4762a1bSJed Brown   */
321c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
322c4762a1bSJed Brown     y = j*hy;
323c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
324c4762a1bSJed Brown       x = i*hx;
32566baab88SBarry Smith       if (PetscApproximateGTE(x,1.0) && PetscApproximateLTE(x,1.5) && PetscApproximateGTE(y,1.0) && PetscApproximateLTE(y,1.5)) u[j][i].v = PetscPowReal(PetscSinReal(4.0*PETSC_PI*x),2.0)*PetscPowReal(PetscSinReal(4.0*PETSC_PI*y),2.0)/4.0;
326c4762a1bSJed Brown       else u[j][i].v = 0.0;
327c4762a1bSJed Brown 
328c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0*u[j][i].v;
329c4762a1bSJed Brown     }
330c4762a1bSJed Brown   }
331c4762a1bSJed Brown 
332c4762a1bSJed Brown   /*
333c4762a1bSJed Brown      Restore vectors
334c4762a1bSJed Brown   */
3359566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,U,&u));
336c4762a1bSJed Brown   PetscFunctionReturn(0);
337c4762a1bSJed Brown }
338c4762a1bSJed Brown 
339c4762a1bSJed Brown PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y)
340c4762a1bSJed Brown {
341c4762a1bSJed Brown    PetscInt i,j,Mx,My,xs,ys,xm,ym;
342c4762a1bSJed Brown    Field **l;
343c4762a1bSJed Brown    PetscFunctionBegin;
344c4762a1bSJed Brown 
3459566063dSJacob 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));
346c4762a1bSJed Brown    /* locate the global i index for x and j index for y */
347c4762a1bSJed Brown    i = (PetscInt)(x*(Mx-1));
348c4762a1bSJed Brown    j = (PetscInt)(y*(My-1));
3499566063dSJacob Faibussowitsch    PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
350c4762a1bSJed Brown 
351c4762a1bSJed Brown    if (xs <= i && i < xs+xm && ys <= j && j < ys+ym) {
352c4762a1bSJed Brown      /* the i,j vertex is on this process */
3539566063dSJacob Faibussowitsch      PetscCall(DMDAVecGetArray(da,lambda,&l));
354c4762a1bSJed Brown      l[j][i].u = 1.0;
355c4762a1bSJed Brown      l[j][i].v = 1.0;
3569566063dSJacob Faibussowitsch      PetscCall(DMDAVecRestoreArray(da,lambda,&l));
357c4762a1bSJed Brown    }
358c4762a1bSJed Brown    PetscFunctionReturn(0);
359c4762a1bSJed Brown }
360c4762a1bSJed Brown 
361c4762a1bSJed Brown PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info,PetscReal t,Field**u,Field**udot,Field**f,void *ptr)
362c4762a1bSJed Brown {
363c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
364c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym;
365c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
366c4762a1bSJed Brown   PetscScalar    uc,uxx,uyy,vc,vxx,vyy;
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   PetscFunctionBegin;
369c4762a1bSJed Brown   hx = 2.50/(PetscReal)(info->mx); sx = 1.0/(hx*hx);
370c4762a1bSJed Brown   hy = 2.50/(PetscReal)(info->my); sy = 1.0/(hy*hy);
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   /* Get local grid boundaries */
373c4762a1bSJed Brown   xs = info->xs; xm = info->xm; ys = info->ys; ym = info->ym;
374c4762a1bSJed Brown 
375c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
376c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
377c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
378c4762a1bSJed Brown       uc        = u[j][i].u;
379c4762a1bSJed Brown       uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
380c4762a1bSJed Brown       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
381c4762a1bSJed Brown       vc        = u[j][i].v;
382c4762a1bSJed Brown       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
383c4762a1bSJed Brown       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
384c4762a1bSJed Brown       f[j][i].u = udot[j][i].u - appctx->D1*(uxx + uyy) + uc*vc*vc - appctx->gamma*(1.0 - uc);
385c4762a1bSJed Brown       f[j][i].v = udot[j][i].v - appctx->D2*(vxx + vyy) - uc*vc*vc + (appctx->gamma + appctx->kappa)*vc;
386c4762a1bSJed Brown     }
387c4762a1bSJed Brown   }
3889566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(16.0*xm*ym));
389c4762a1bSJed Brown   PetscFunctionReturn(0);
390c4762a1bSJed Brown }
391c4762a1bSJed Brown 
392c4762a1bSJed Brown PetscErrorCode IFunctionActive(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
393c4762a1bSJed Brown {
394c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
395c4762a1bSJed Brown   DM             da;
396c4762a1bSJed Brown   DMDALocalInfo  info;
397c4762a1bSJed Brown   Field          **u,**f,**udot;
398c4762a1bSJed Brown   Vec            localU;
399c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,gxs,gys,gxm,gym;
400c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
401c4762a1bSJed Brown   adouble        uc,uxx,uyy,vc,vxx,vyy;
402c4762a1bSJed Brown   AField         **f_a,*f_c,**u_a,*u_c;
403c4762a1bSJed Brown   PetscScalar    dummy;
404c4762a1bSJed Brown 
405c4762a1bSJed Brown   PetscFunctionBegin;
406c4762a1bSJed Brown 
4079566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
4089566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
4099566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
410c4762a1bSJed Brown   hx = 2.50/(PetscReal)(info.mx); sx = 1.0/(hx*hx);
411c4762a1bSJed Brown   hy = 2.50/(PetscReal)(info.my); sy = 1.0/(hy*hy);
412c4762a1bSJed Brown   xs = info.xs; xm = info.xm; gxs = info.gxs; gxm = info.gxm;
413c4762a1bSJed Brown   ys = info.ys; ym = info.ym; gys = info.gys; gym = info.gym;
414c4762a1bSJed Brown 
415c4762a1bSJed Brown   /*
416c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
417c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
418c4762a1bSJed Brown      By placing code between these two statements, computations can be
419c4762a1bSJed Brown      done while messages are in transition.
420c4762a1bSJed Brown   */
4219566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
4229566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   /*
425c4762a1bSJed Brown      Get pointers to vector data
426c4762a1bSJed Brown   */
4279566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
4289566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,F,&f));
4299566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,Udot,&udot));
430c4762a1bSJed Brown 
431c4762a1bSJed Brown   /*
432c4762a1bSJed Brown     Create contiguous 1-arrays of AFields
433c4762a1bSJed Brown 
434c4762a1bSJed Brown     NOTE: Memory for ADOL-C active variables (such as adouble and AField)
435c4762a1bSJed Brown           cannot be allocated using PetscMalloc, as this does not call the
436c4762a1bSJed Brown           relevant class constructor. Instead, we use the C++ keyword `new`.
437c4762a1bSJed Brown   */
438c4762a1bSJed Brown   u_c = new AField[info.gxm*info.gym];
439c4762a1bSJed Brown   f_c = new AField[info.gxm*info.gym];
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   /* Create corresponding 2-arrays of AFields */
442c4762a1bSJed Brown   u_a = new AField*[info.gym];
443c4762a1bSJed Brown   f_a = new AField*[info.gym];
444c4762a1bSJed Brown 
445c4762a1bSJed Brown   /* Align indices between array types to endow 2d array with ghost points */
4469566063dSJacob Faibussowitsch   PetscCall(GiveGhostPoints(da,u_c,&u_a));
4479566063dSJacob Faibussowitsch   PetscCall(GiveGhostPoints(da,f_c,&f_a));
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   trace_on(1);  /* Start of active section on tape 1 */
450c4762a1bSJed Brown 
451c4762a1bSJed Brown   /*
452c4762a1bSJed Brown     Mark independence
453c4762a1bSJed Brown 
454c4762a1bSJed Brown     NOTE: Ghost points are marked as independent, in place of the points they represent on
455c4762a1bSJed Brown           other processors / on other boundaries.
456c4762a1bSJed Brown   */
457c4762a1bSJed Brown   for (j=gys; j<gys+gym; j++) {
458c4762a1bSJed Brown     for (i=gxs; i<gxs+gxm; i++) {
459c4762a1bSJed Brown       u_a[j][i].u <<= u[j][i].u;
460c4762a1bSJed Brown       u_a[j][i].v <<= u[j][i].v;
461c4762a1bSJed Brown     }
462c4762a1bSJed Brown   }
463c4762a1bSJed Brown 
464c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
465c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
466c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
467c4762a1bSJed Brown       uc        = u_a[j][i].u;
468c4762a1bSJed Brown       uxx       = (-2.0*uc + u_a[j][i-1].u + u_a[j][i+1].u)*sx;
469c4762a1bSJed Brown       uyy       = (-2.0*uc + u_a[j-1][i].u + u_a[j+1][i].u)*sy;
470c4762a1bSJed Brown       vc        = u_a[j][i].v;
471c4762a1bSJed Brown       vxx       = (-2.0*vc + u_a[j][i-1].v + u_a[j][i+1].v)*sx;
472c4762a1bSJed Brown       vyy       = (-2.0*vc + u_a[j-1][i].v + u_a[j+1][i].v)*sy;
473c4762a1bSJed Brown       f_a[j][i].u = udot[j][i].u - appctx->D1*(uxx + uyy) + uc*vc*vc - appctx->gamma*(1.0 - uc);
474c4762a1bSJed Brown       f_a[j][i].v = udot[j][i].v - appctx->D2*(vxx + vyy) - uc*vc*vc + (appctx->gamma + appctx->kappa)*vc;
475c4762a1bSJed Brown     }
476c4762a1bSJed Brown   }
477c4762a1bSJed Brown 
478c4762a1bSJed Brown   /*
479c4762a1bSJed Brown     Mark dependence
480c4762a1bSJed Brown 
481c4762a1bSJed Brown     NOTE: Marking dependence of dummy variables makes the index notation much simpler when forming
482c4762a1bSJed Brown           the Jacobian later.
483c4762a1bSJed Brown   */
484c4762a1bSJed Brown   for (j=gys; j<gys+gym; j++) {
485c4762a1bSJed Brown     for (i=gxs; i<gxs+gxm; i++) {
486c4762a1bSJed Brown       if ((i < xs) || (i >= xs+xm) || (j < ys) || (j >= ys+ym)) {
487c4762a1bSJed Brown         f_a[j][i].u >>= dummy;
488c4762a1bSJed Brown         f_a[j][i].v >>= dummy;
489c4762a1bSJed Brown       } else {
490c4762a1bSJed Brown         f_a[j][i].u >>= f[j][i].u;
491c4762a1bSJed Brown         f_a[j][i].v >>= f[j][i].v;
492c4762a1bSJed Brown       }
493c4762a1bSJed Brown     }
494c4762a1bSJed Brown   }
495c4762a1bSJed Brown   trace_off();  /* End of active section */
4969566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(16.0*xm*ym));
497c4762a1bSJed Brown 
498c4762a1bSJed Brown   /* Restore vectors */
4999566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,F,&f));
5009566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
5019566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,Udot,&udot));
502c4762a1bSJed Brown 
5039566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
504410585f6SHong Zhang 
505c4762a1bSJed Brown   /* Destroy AFields appropriately */
506c4762a1bSJed Brown   f_a += info.gys;
507c4762a1bSJed Brown   u_a += info.gys;
508c4762a1bSJed Brown   delete[] f_a;
509c4762a1bSJed Brown   delete[] u_a;
510c4762a1bSJed Brown   delete[] f_c;
511c4762a1bSJed Brown   delete[] u_c;
512c4762a1bSJed Brown 
513c4762a1bSJed Brown   PetscFunctionReturn(0);
514c4762a1bSJed Brown }
515c4762a1bSJed Brown 
516c4762a1bSJed Brown PetscErrorCode RHSFunctionPassive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr)
517c4762a1bSJed Brown {
518c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
519c4762a1bSJed Brown   DM             da;
520c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
521c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
522c4762a1bSJed Brown   PetscScalar    uc,uxx,uyy,vc,vxx,vyy;
523c4762a1bSJed Brown   Field          **u,**f;
524c4762a1bSJed Brown   Vec            localU,localF;
525c4762a1bSJed Brown 
526c4762a1bSJed Brown   PetscFunctionBegin;
5279566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
5289566063dSJacob 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));
529c4762a1bSJed Brown   hx = 2.50/(PetscReal)(Mx);sx = 1.0/(hx*hx);
530c4762a1bSJed Brown   hy = 2.50/(PetscReal)(My);sy = 1.0/(hy*hy);
5319566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
5329566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localF));
533c4762a1bSJed Brown 
534c4762a1bSJed Brown   /*
535c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
536c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
537c4762a1bSJed Brown      By placing code between these two statements, computations can be
538c4762a1bSJed Brown      done while messages are in transition.
539c4762a1bSJed Brown   */
5409566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
5419566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
5429566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F)); // NOTE (1): See (2) below
5439566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,F,INSERT_VALUES,localF));
5449566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,F,INSERT_VALUES,localF));
545c4762a1bSJed Brown 
546c4762a1bSJed Brown   /*
547c4762a1bSJed Brown      Get pointers to vector data
548c4762a1bSJed Brown   */
5499566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
5509566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,localF,&f));
551c4762a1bSJed Brown 
552c4762a1bSJed Brown   /*
553c4762a1bSJed Brown      Get local grid boundaries
554c4762a1bSJed Brown   */
5559566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
556c4762a1bSJed Brown 
557c4762a1bSJed Brown   /*
558c4762a1bSJed Brown      Compute function over the locally owned part of the grid
559c4762a1bSJed Brown   */
560c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
561c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
562c4762a1bSJed Brown       uc        = u[j][i].u;
563c4762a1bSJed Brown       uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
564c4762a1bSJed Brown       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
565c4762a1bSJed Brown       vc        = u[j][i].v;
566c4762a1bSJed Brown       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
567c4762a1bSJed Brown       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
568c4762a1bSJed Brown       f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);
569c4762a1bSJed Brown       f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc;
570c4762a1bSJed Brown     }
571c4762a1bSJed Brown   }
572c4762a1bSJed Brown 
573c4762a1bSJed Brown   /*
574c4762a1bSJed Brown      Gather global vector, using the 2-step process
575c4762a1bSJed Brown         DMLocalToGlobalBegin(),DMLocalToGlobalEnd().
576c4762a1bSJed Brown 
577c4762a1bSJed Brown      NOTE (2): We need to use ADD_VALUES if boundaries are not of type DM_BOUNDARY_NONE or
578c4762a1bSJed Brown                DM_BOUNDARY_GHOSTED, meaning we should also zero F before addition (see (1) above).
579c4762a1bSJed Brown   */
5809566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da,localF,ADD_VALUES,F));
5819566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da,localF,ADD_VALUES,F));
582c4762a1bSJed Brown 
583c4762a1bSJed Brown   /*
584c4762a1bSJed Brown      Restore vectors
585c4762a1bSJed Brown   */
5869566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,localF,&f));
5879566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
5889566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localF));
5899566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
5909566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(16.0*xm*ym));
591c4762a1bSJed Brown   PetscFunctionReturn(0);
592c4762a1bSJed Brown }
593c4762a1bSJed Brown 
594c4762a1bSJed Brown PetscErrorCode RHSFunctionActive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr)
595c4762a1bSJed Brown {
596c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
597c4762a1bSJed Brown   DM             da;
598c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,gxs,gys,gxm,gym,Mx,My;
599c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
600c4762a1bSJed Brown   AField         **f_a,*f_c,**u_a,*u_c;
601c4762a1bSJed Brown   adouble        uc,uxx,uyy,vc,vxx,vyy;
602c4762a1bSJed Brown   Field          **u,**f;
603c4762a1bSJed Brown   Vec            localU,localF;
604c4762a1bSJed Brown 
605c4762a1bSJed Brown   PetscFunctionBegin;
6069566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
6079566063dSJacob 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));
608c4762a1bSJed Brown   hx = 2.50/(PetscReal)(Mx);sx = 1.0/(hx*hx);
609c4762a1bSJed Brown   hy = 2.50/(PetscReal)(My);sy = 1.0/(hy*hy);
6109566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
6119566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localF));
612c4762a1bSJed Brown 
613c4762a1bSJed Brown   /*
614c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
615c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
616c4762a1bSJed Brown      By placing code between these two statements, computations can be
617c4762a1bSJed Brown      done while messages are in transition.
618c4762a1bSJed Brown   */
6199566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
6209566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
6219566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F)); // NOTE (1): See (2) below
6229566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,F,INSERT_VALUES,localF));
6239566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,F,INSERT_VALUES,localF));
624c4762a1bSJed Brown 
625c4762a1bSJed Brown   /*
626c4762a1bSJed Brown      Get pointers to vector data
627c4762a1bSJed Brown   */
6289566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
6299566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,localF,&f));
630c4762a1bSJed Brown 
631c4762a1bSJed Brown   /*
632c4762a1bSJed Brown      Get local and ghosted grid boundaries
633c4762a1bSJed Brown   */
6349566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,NULL,&gxm,&gym,NULL));
6359566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
636c4762a1bSJed Brown 
637c4762a1bSJed Brown   /*
638c4762a1bSJed Brown     Create contiguous 1-arrays of AFields
639c4762a1bSJed Brown 
640c4762a1bSJed Brown     NOTE: Memory for ADOL-C active variables (such as adouble and AField)
641c4762a1bSJed Brown           cannot be allocated using PetscMalloc, as this does not call the
642c4762a1bSJed Brown           relevant class constructor. Instead, we use the C++ keyword `new`.
643c4762a1bSJed Brown   */
644c4762a1bSJed Brown   u_c = new AField[gxm*gym];
645c4762a1bSJed Brown   f_c = new AField[gxm*gym];
646c4762a1bSJed Brown 
647c4762a1bSJed Brown   /* Create corresponding 2-arrays of AFields */
648c4762a1bSJed Brown   u_a = new AField*[gym];
649c4762a1bSJed Brown   f_a = new AField*[gym];
650c4762a1bSJed Brown 
651c4762a1bSJed Brown   /* Align indices between array types to endow 2d array with ghost points */
6529566063dSJacob Faibussowitsch   PetscCall(GiveGhostPoints(da,u_c,&u_a));
6539566063dSJacob Faibussowitsch   PetscCall(GiveGhostPoints(da,f_c,&f_a));
654c4762a1bSJed Brown 
655c4762a1bSJed Brown   /*
656c4762a1bSJed Brown      Compute function over the locally owned part of the grid
657c4762a1bSJed Brown   */
658c4762a1bSJed Brown   trace_on(1);  // ----------------------------------------------- Start of active section
659c4762a1bSJed Brown 
660c4762a1bSJed Brown   /*
661c4762a1bSJed Brown     Mark independence
662c4762a1bSJed Brown 
663c4762a1bSJed Brown     NOTE: Ghost points are marked as independent, in place of the points they represent on
664c4762a1bSJed Brown           other processors / on other boundaries.
665c4762a1bSJed Brown   */
666c4762a1bSJed Brown   for (j=gys; j<gys+gym; j++) {
667c4762a1bSJed Brown     for (i=gxs; i<gxs+gxm; i++) {
668c4762a1bSJed Brown       u_a[j][i].u <<= u[j][i].u;
669c4762a1bSJed Brown       u_a[j][i].v <<= u[j][i].v;
670c4762a1bSJed Brown     }
671c4762a1bSJed Brown   }
672c4762a1bSJed Brown 
673c4762a1bSJed Brown   /*
674c4762a1bSJed Brown     Compute function over the locally owned part of the grid
675c4762a1bSJed Brown   */
676c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
677c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
678c4762a1bSJed Brown       uc          = u_a[j][i].u;
679c4762a1bSJed Brown       uxx         = (-2.0*uc + u_a[j][i-1].u + u_a[j][i+1].u)*sx;
680c4762a1bSJed Brown       uyy         = (-2.0*uc + u_a[j-1][i].u + u_a[j+1][i].u)*sy;
681c4762a1bSJed Brown       vc          = u_a[j][i].v;
682c4762a1bSJed Brown       vxx         = (-2.0*vc + u_a[j][i-1].v + u_a[j][i+1].v)*sx;
683c4762a1bSJed Brown       vyy         = (-2.0*vc + u_a[j-1][i].v + u_a[j+1][i].v)*sy;
684c4762a1bSJed Brown       f_a[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);
685c4762a1bSJed Brown       f_a[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc;
686c4762a1bSJed Brown     }
687c4762a1bSJed Brown   }
688c4762a1bSJed Brown   /*
689c4762a1bSJed Brown     Mark dependence
690c4762a1bSJed Brown 
691c4762a1bSJed Brown     NOTE: Ghost points are marked as dependent in order to vastly simplify index notation
692c4762a1bSJed Brown           during Jacobian assembly.
693c4762a1bSJed Brown   */
694c4762a1bSJed Brown   for (j=gys; j<gys+gym; j++) {
695c4762a1bSJed Brown     for (i=gxs; i<gxs+gxm; i++) {
696c4762a1bSJed Brown       f_a[j][i].u >>= f[j][i].u;
697c4762a1bSJed Brown       f_a[j][i].v >>= f[j][i].v;
698c4762a1bSJed Brown     }
699c4762a1bSJed Brown   }
700c4762a1bSJed Brown   trace_off();  // ----------------------------------------------- End of active section
701c4762a1bSJed Brown 
702c4762a1bSJed Brown   /* Test zeroth order scalar evaluation in ADOL-C gives the same result */
703c4762a1bSJed Brown //  if (appctx->adctx->zos) {
7049566063dSJacob Faibussowitsch //    PetscCall(TestZOS2d(da,f,u,appctx)); // FIXME: Why does this give nonzero?
705c4762a1bSJed Brown //  }
706c4762a1bSJed Brown 
707c4762a1bSJed Brown   /*
708c4762a1bSJed Brown      Gather global vector, using the 2-step process
709c4762a1bSJed Brown         DMLocalToGlobalBegin(),DMLocalToGlobalEnd().
710c4762a1bSJed Brown 
711c4762a1bSJed Brown      NOTE (2): We need to use ADD_VALUES if boundaries are not of type DM_BOUNDARY_NONE or
712c4762a1bSJed Brown                DM_BOUNDARY_GHOSTED, meaning we should also zero F before addition (see (1) above).
713c4762a1bSJed Brown   */
7149566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(da,localF,ADD_VALUES,F));
7159566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(da,localF,ADD_VALUES,F));
716c4762a1bSJed Brown 
717c4762a1bSJed Brown   /*
718c4762a1bSJed Brown      Restore vectors
719c4762a1bSJed Brown   */
7209566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,localF,&f));
7219566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
7229566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localF));
7239566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
7249566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(16.0*xm*ym));
725c4762a1bSJed Brown 
726c4762a1bSJed Brown   /* Destroy AFields appropriately */
727c4762a1bSJed Brown   f_a += gys;
728c4762a1bSJed Brown   u_a += gys;
729c4762a1bSJed Brown   delete[] f_a;
730c4762a1bSJed Brown   delete[] u_a;
731c4762a1bSJed Brown   delete[] f_c;
732c4762a1bSJed Brown   delete[] u_c;
733c4762a1bSJed Brown 
734c4762a1bSJed Brown   PetscFunctionReturn(0);
735c4762a1bSJed Brown }
736c4762a1bSJed Brown 
737c4762a1bSJed Brown PetscErrorCode IJacobianAdolc(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
738c4762a1bSJed Brown {
739c4762a1bSJed Brown   AppCtx            *appctx = (AppCtx*)ctx;
740c4762a1bSJed Brown   DM                da;
741a8c08197SHong Zhang   const PetscScalar *u_vec;
742c4762a1bSJed Brown   Vec               localU;
743c4762a1bSJed Brown 
744c4762a1bSJed Brown   PetscFunctionBegin;
7459566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
7469566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
747c4762a1bSJed Brown 
748c4762a1bSJed Brown   /*
749c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
750c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
751c4762a1bSJed Brown      By placing code between these two statements, computations can be
752c4762a1bSJed Brown      done while messages are in transition.
753c4762a1bSJed Brown   */
7549566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
7559566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
756c4762a1bSJed Brown 
757c4762a1bSJed Brown   /* Get pointers to vector data */
7589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(localU,&u_vec));
759c4762a1bSJed Brown 
760c4762a1bSJed Brown   /*
761c4762a1bSJed Brown     Compute Jacobian
762c4762a1bSJed Brown   */
7639566063dSJacob Faibussowitsch   PetscCall(PetscAdolcComputeIJacobianLocalIDMass(1,A,u_vec,a,appctx->adctx));
764c4762a1bSJed Brown 
765c4762a1bSJed Brown   /*
766c4762a1bSJed Brown      Restore vectors
767c4762a1bSJed Brown   */
7689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(localU,&u_vec));
7699566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
770c4762a1bSJed Brown   PetscFunctionReturn(0);
771c4762a1bSJed Brown }
772c4762a1bSJed Brown 
773c4762a1bSJed Brown PetscErrorCode IJacobianByHand(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
774c4762a1bSJed Brown {
775c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
776c4762a1bSJed Brown   DM             da;
777c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
778c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
779c4762a1bSJed Brown   PetscScalar    uc,vc;
780c4762a1bSJed Brown   Field          **u;
781c4762a1bSJed Brown   Vec            localU;
782c4762a1bSJed Brown   MatStencil     stencil[6],rowstencil;
783c4762a1bSJed Brown   PetscScalar    entries[6];
784c4762a1bSJed Brown 
785c4762a1bSJed Brown   PetscFunctionBegin;
7869566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
7879566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
7889566063dSJacob 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));
789c4762a1bSJed Brown 
790c4762a1bSJed Brown   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
791c4762a1bSJed Brown   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
792c4762a1bSJed Brown 
793c4762a1bSJed Brown   /*
794c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
795c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
796c4762a1bSJed Brown      By placing code between these two statements, computations can be
797c4762a1bSJed Brown      done while messages are in transition.
798c4762a1bSJed Brown   */
7999566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
8009566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
801c4762a1bSJed Brown 
802c4762a1bSJed Brown   /*
803c4762a1bSJed Brown      Get pointers to vector data
804c4762a1bSJed Brown   */
8059566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
806c4762a1bSJed Brown 
807c4762a1bSJed Brown   /*
808c4762a1bSJed Brown      Get local grid boundaries
809c4762a1bSJed Brown   */
8109566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
811c4762a1bSJed Brown 
812c4762a1bSJed Brown   stencil[0].k = 0;
813c4762a1bSJed Brown   stencil[1].k = 0;
814c4762a1bSJed Brown   stencil[2].k = 0;
815c4762a1bSJed Brown   stencil[3].k = 0;
816c4762a1bSJed Brown   stencil[4].k = 0;
817c4762a1bSJed Brown   stencil[5].k = 0;
818c4762a1bSJed Brown   rowstencil.k = 0;
819c4762a1bSJed Brown   /*
820c4762a1bSJed Brown      Compute function over the locally owned part of the grid
821c4762a1bSJed Brown   */
822c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
823c4762a1bSJed Brown 
824c4762a1bSJed Brown     stencil[0].j = j-1;
825c4762a1bSJed Brown     stencil[1].j = j+1;
826c4762a1bSJed Brown     stencil[2].j = j;
827c4762a1bSJed Brown     stencil[3].j = j;
828c4762a1bSJed Brown     stencil[4].j = j;
829c4762a1bSJed Brown     stencil[5].j = j;
830c4762a1bSJed Brown     rowstencil.k = 0; rowstencil.j = j;
831c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
832c4762a1bSJed Brown       uc = u[j][i].u;
833c4762a1bSJed Brown       vc = u[j][i].v;
834c4762a1bSJed Brown 
835c4762a1bSJed Brown       /*      uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
836c4762a1bSJed Brown       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
837c4762a1bSJed Brown 
838c4762a1bSJed Brown       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
839c4762a1bSJed Brown       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
840c4762a1bSJed Brown        f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
841c4762a1bSJed Brown 
842c4762a1bSJed Brown       stencil[0].i = i; stencil[0].c = 0; entries[0] = -appctx->D1*sy;
843c4762a1bSJed Brown       stencil[1].i = i; stencil[1].c = 0; entries[1] = -appctx->D1*sy;
844c4762a1bSJed Brown       stencil[2].i = i-1; stencil[2].c = 0; entries[2] = -appctx->D1*sx;
845c4762a1bSJed Brown       stencil[3].i = i+1; stencil[3].c = 0; entries[3] = -appctx->D1*sx;
846c4762a1bSJed Brown       stencil[4].i = i; stencil[4].c = 0; entries[4] = 2.0*appctx->D1*(sx + sy) + vc*vc + appctx->gamma + a;
847c4762a1bSJed Brown       stencil[5].i = i; stencil[5].c = 1; entries[5] = 2.0*uc*vc;
848c4762a1bSJed Brown       rowstencil.i = i; rowstencil.c = 0;
849c4762a1bSJed Brown 
8509566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
851c4762a1bSJed Brown       if (appctx->aijpc) {
8529566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
853c4762a1bSJed Brown       }
854c4762a1bSJed Brown       stencil[0].c = 1; entries[0] = -appctx->D2*sy;
855c4762a1bSJed Brown       stencil[1].c = 1; entries[1] = -appctx->D2*sy;
856c4762a1bSJed Brown       stencil[2].c = 1; entries[2] = -appctx->D2*sx;
857c4762a1bSJed Brown       stencil[3].c = 1; entries[3] = -appctx->D2*sx;
858c4762a1bSJed Brown       stencil[4].c = 1; entries[4] = 2.0*appctx->D2*(sx + sy) - 2.0*uc*vc + appctx->gamma + appctx->kappa + a;
859c4762a1bSJed Brown       stencil[5].c = 0; entries[5] = -vc*vc;
860c4762a1bSJed Brown 
8619566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
862c4762a1bSJed Brown       if (appctx->aijpc) {
8639566063dSJacob Faibussowitsch         PetscCall(MatSetValuesStencil(B,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
864c4762a1bSJed Brown       }
865c4762a1bSJed Brown       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
866c4762a1bSJed Brown     }
867c4762a1bSJed Brown   }
868c4762a1bSJed Brown 
869c4762a1bSJed Brown   /*
870c4762a1bSJed Brown      Restore vectors
871c4762a1bSJed Brown   */
8729566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(19.0*xm*ym));
8739566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
8749566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
8759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
8769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
8779566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
878c4762a1bSJed Brown   if (appctx->aijpc) {
8799566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
8809566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
8819566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
882c4762a1bSJed Brown   }
883c4762a1bSJed Brown   PetscFunctionReturn(0);
884c4762a1bSJed Brown }
885c4762a1bSJed Brown 
886c4762a1bSJed Brown PetscErrorCode RHSJacobianByHand(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
887c4762a1bSJed Brown {
888c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
889c4762a1bSJed Brown   DM             da;
890c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
891c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
892c4762a1bSJed Brown   PetscScalar    uc,vc;
893c4762a1bSJed Brown   Field          **u;
894c4762a1bSJed Brown   Vec            localU;
895c4762a1bSJed Brown   MatStencil     stencil[6],rowstencil;
896c4762a1bSJed Brown   PetscScalar    entries[6];
897c4762a1bSJed Brown 
898c4762a1bSJed Brown   PetscFunctionBegin;
8999566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
9009566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
9019566063dSJacob 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));
902c4762a1bSJed Brown 
903c4762a1bSJed Brown   hx = 2.50/(PetscReal)(Mx); sx = 1.0/(hx*hx);
904c4762a1bSJed Brown   hy = 2.50/(PetscReal)(My); sy = 1.0/(hy*hy);
905c4762a1bSJed Brown 
906c4762a1bSJed Brown   /*
907c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
908c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
909c4762a1bSJed Brown      By placing code between these two statements, computations can be
910c4762a1bSJed Brown      done while messages are in transition.
911c4762a1bSJed Brown   */
9129566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
9139566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
914c4762a1bSJed Brown 
915c4762a1bSJed Brown   /*
916c4762a1bSJed Brown      Get pointers to vector data
917c4762a1bSJed Brown   */
9189566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
919c4762a1bSJed Brown 
920c4762a1bSJed Brown   /*
921c4762a1bSJed Brown      Get local grid boundaries
922c4762a1bSJed Brown   */
9239566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
924c4762a1bSJed Brown 
925c4762a1bSJed Brown   stencil[0].k = 0;
926c4762a1bSJed Brown   stencil[1].k = 0;
927c4762a1bSJed Brown   stencil[2].k = 0;
928c4762a1bSJed Brown   stencil[3].k = 0;
929c4762a1bSJed Brown   stencil[4].k = 0;
930c4762a1bSJed Brown   stencil[5].k = 0;
931c4762a1bSJed Brown   rowstencil.k = 0;
932c4762a1bSJed Brown 
933c4762a1bSJed Brown   /*
934c4762a1bSJed Brown      Compute function over the locally owned part of the grid
935c4762a1bSJed Brown   */
936c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
937c4762a1bSJed Brown 
938c4762a1bSJed Brown     stencil[0].j = j-1;
939c4762a1bSJed Brown     stencil[1].j = j+1;
940c4762a1bSJed Brown     stencil[2].j = j;
941c4762a1bSJed Brown     stencil[3].j = j;
942c4762a1bSJed Brown     stencil[4].j = j;
943c4762a1bSJed Brown     stencil[5].j = j;
944c4762a1bSJed Brown     rowstencil.k = 0; rowstencil.j = j;
945c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
946c4762a1bSJed Brown       uc = u[j][i].u;
947c4762a1bSJed Brown       vc = u[j][i].v;
948c4762a1bSJed Brown 
949c4762a1bSJed Brown       /*      uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
950c4762a1bSJed Brown       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
951c4762a1bSJed Brown 
952c4762a1bSJed Brown       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
953c4762a1bSJed Brown       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
954c4762a1bSJed Brown        f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
955c4762a1bSJed Brown 
956c4762a1bSJed Brown       stencil[0].i = i; stencil[0].c = 0; entries[0] = appctx->D1*sy;
957c4762a1bSJed Brown       stencil[1].i = i; stencil[1].c = 0; entries[1] = appctx->D1*sy;
958c4762a1bSJed Brown       stencil[2].i = i-1; stencil[2].c = 0; entries[2] = appctx->D1*sx;
959c4762a1bSJed Brown       stencil[3].i = i+1; stencil[3].c = 0; entries[3] = appctx->D1*sx;
960c4762a1bSJed Brown       stencil[4].i = i; stencil[4].c = 0; entries[4] = -2.0*appctx->D1*(sx + sy) - vc*vc - appctx->gamma;
961c4762a1bSJed Brown       stencil[5].i = i; stencil[5].c = 1; entries[5] = -2.0*uc*vc;
962c4762a1bSJed Brown       rowstencil.i = i; rowstencil.c = 0;
963c4762a1bSJed Brown 
9649566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
965c4762a1bSJed Brown 
966c4762a1bSJed Brown       stencil[0].c = 1; entries[0] = appctx->D2*sy;
967c4762a1bSJed Brown       stencil[1].c = 1; entries[1] = appctx->D2*sy;
968c4762a1bSJed Brown       stencil[2].c = 1; entries[2] = appctx->D2*sx;
969c4762a1bSJed Brown       stencil[3].c = 1; entries[3] = appctx->D2*sx;
970c4762a1bSJed Brown       stencil[4].c = 1; entries[4] = -2.0*appctx->D2*(sx + sy) + 2.0*uc*vc - appctx->gamma - appctx->kappa;
971c4762a1bSJed Brown       stencil[5].c = 0; entries[5] = vc*vc;
972c4762a1bSJed Brown       rowstencil.c = 1;
973c4762a1bSJed Brown 
9749566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
975c4762a1bSJed Brown       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
976c4762a1bSJed Brown     }
977c4762a1bSJed Brown   }
978c4762a1bSJed Brown 
979c4762a1bSJed Brown   /*
980c4762a1bSJed Brown      Restore vectors
981c4762a1bSJed Brown   */
9829566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(19.0*xm*ym));
9839566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
9849566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
9859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
9869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
9879566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
988c4762a1bSJed Brown   if (appctx->aijpc) {
9899566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
9909566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
9919566063dSJacob Faibussowitsch     PetscCall(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
992c4762a1bSJed Brown   }
993c4762a1bSJed Brown   PetscFunctionReturn(0);
994c4762a1bSJed Brown }
995c4762a1bSJed Brown 
996c4762a1bSJed Brown PetscErrorCode RHSJacobianAdolc(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
997c4762a1bSJed Brown {
998c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;
999c4762a1bSJed Brown   DM             da;
1000c4762a1bSJed Brown   PetscScalar    *u_vec;
1001c4762a1bSJed Brown   Vec            localU;
1002c4762a1bSJed Brown 
1003c4762a1bSJed Brown   PetscFunctionBegin;
10049566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
10059566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
1006c4762a1bSJed Brown 
1007c4762a1bSJed Brown   /*
1008c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
1009c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
1010c4762a1bSJed Brown      By placing code between these two statements, computations can be
1011c4762a1bSJed Brown      done while messages are in transition.
1012c4762a1bSJed Brown   */
10139566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
10149566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
1015c4762a1bSJed Brown 
1016c4762a1bSJed Brown   /* Get pointers to vector data */
10179566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localU,&u_vec));
1018c4762a1bSJed Brown 
1019c4762a1bSJed Brown   /*
1020c4762a1bSJed Brown     Compute Jacobian
1021c4762a1bSJed Brown   */
10229566063dSJacob Faibussowitsch   PetscCall(PetscAdolcComputeRHSJacobianLocal(1,A,u_vec,appctx->adctx));
1023c4762a1bSJed Brown 
1024c4762a1bSJed Brown   /*
1025c4762a1bSJed Brown      Restore vectors
1026c4762a1bSJed Brown   */
10279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localU,&u_vec));
10289566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
1029c4762a1bSJed Brown   PetscFunctionReturn(0);
1030c4762a1bSJed Brown }
1031c4762a1bSJed Brown 
1032c4762a1bSJed Brown /*TEST
1033c4762a1bSJed Brown 
1034c4762a1bSJed Brown   build:
1035c4762a1bSJed Brown     requires: double !complex adolc colpack
1036c4762a1bSJed Brown 
1037c4762a1bSJed Brown   test:
1038c4762a1bSJed Brown     suffix: 1
1039c4762a1bSJed Brown     nsize: 1
1040c4762a1bSJed Brown     args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian
1041c4762a1bSJed Brown     output_file: output/adr_ex5adj_1.out
1042c4762a1bSJed Brown 
1043c4762a1bSJed Brown   test:
1044c4762a1bSJed Brown     suffix: 2
1045c4762a1bSJed Brown     nsize: 1
1046c4762a1bSJed Brown     args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian -implicitform
1047c4762a1bSJed Brown     output_file: output/adr_ex5adj_2.out
1048c4762a1bSJed Brown 
1049c4762a1bSJed Brown   test:
1050c4762a1bSJed Brown     suffix: 3
1051c4762a1bSJed Brown     nsize: 4
1052c4762a1bSJed Brown     args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor
1053c4762a1bSJed Brown     output_file: output/adr_ex5adj_3.out
1054c4762a1bSJed Brown 
1055c4762a1bSJed Brown   test:
1056c4762a1bSJed Brown     suffix: 4
1057c4762a1bSJed Brown     nsize: 4
1058c4762a1bSJed Brown     args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor -implicitform
1059c4762a1bSJed Brown     output_file: output/adr_ex5adj_4.out
1060c4762a1bSJed Brown 
1061c4762a1bSJed Brown   testset:
1062c4762a1bSJed Brown     suffix: 5
1063c4762a1bSJed Brown     nsize: 4
1064c4762a1bSJed Brown     args: -ts_max_steps 10 -da_grid_x 15 -da_grid_y 15 -ts_monitor -ts_adjoint_monitor -adolc_sparse
1065c4762a1bSJed Brown     output_file: output/adr_ex5adj_5.out
1066c4762a1bSJed Brown 
1067c4762a1bSJed Brown   testset:
1068c4762a1bSJed Brown     suffix: 6
1069c4762a1bSJed Brown     nsize: 4
1070c4762a1bSJed Brown     args: -ts_max_steps 10 -da_grid_x 15 -da_grid_y 15 -ts_monitor -ts_adjoint_monitor -adolc_sparse -implicitform
1071c4762a1bSJed Brown     output_file: output/adr_ex5adj_6.out
1072c4762a1bSJed Brown 
1073c4762a1bSJed Brown TEST*/
1074