xref: /petsc/src/ts/tutorials/autodiff/adr_ex5adj.cxx (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Demonstrates adjoint sensitivity analysis for Reaction-Diffusion Equations.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    Concepts: TS^time-dependent nonlinear problems
5c4762a1bSJed Brown    Concepts: Automatic differentiation using ADOL-C
6c4762a1bSJed Brown */
7c4762a1bSJed Brown /*
8c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-adolc.
9c4762a1bSJed Brown 
10c4762a1bSJed Brown    For documentation on ADOL-C, see
11c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
12c4762a1bSJed Brown 
13c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-colpack
14c4762a1bSJed Brown 
15c4762a1bSJed Brown    For documentation on ColPack, see
16c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/git.colpack/README.md
17c4762a1bSJed Brown */
18c4762a1bSJed Brown /* ------------------------------------------------------------------------
19c4762a1bSJed Brown   See ../advection-diffusion-reaction/ex5 for a description of the problem
20c4762a1bSJed Brown   ------------------------------------------------------------------------- */
21c4762a1bSJed Brown 
22c4762a1bSJed Brown /*
23c4762a1bSJed Brown   Runtime options:
24c4762a1bSJed Brown 
25c4762a1bSJed Brown     Solver:
26c4762a1bSJed Brown       -forwardonly       - Run the forward simulation without adjoint.
27c4762a1bSJed Brown       -implicitform      - Provide IFunction and IJacobian to TS, if not set, RHSFunction and RHSJacobian will be used.
28c4762a1bSJed Brown       -aijpc             - Set the preconditioner matrix to be aij (the Jacobian matrix can be of a different type such as ELL).
29c4762a1bSJed Brown 
30c4762a1bSJed Brown     Jacobian generation:
31c4762a1bSJed Brown       -jacobian_by_hand  - Use the hand-coded Jacobian of ex5.c, rather than generating it automatically.
32c4762a1bSJed Brown       -no_annotation     - Do not annotate ADOL-C active variables. (Should be used alongside -jacobian_by_hand.)
33c4762a1bSJed Brown       -adolc_sparse      - Calculate Jacobian in compressed format, using ADOL-C sparse functionality.
34c4762a1bSJed Brown       -adolc_sparse_view - Print sparsity pattern used by -adolc_sparse option.
35c4762a1bSJed Brown  */
36c4762a1bSJed Brown /*
37c4762a1bSJed Brown   NOTE: If -adolc_sparse option is used, at least four processors should be used, so that no processor owns boundaries which are
38c4762a1bSJed Brown         identified by the periodic boundary conditions. The number of grid points in both x- and y-directions should be multiples
39c4762a1bSJed Brown         of 5, in order for the 5-point stencil to be cleanly parallelised.
40c4762a1bSJed Brown */
41c4762a1bSJed Brown 
42c4762a1bSJed Brown #include <petscdmda.h>
43c4762a1bSJed Brown #include <petscts.h>
44c4762a1bSJed Brown #include "adolc-utils/drivers.cxx"
45c4762a1bSJed Brown #include <adolc/adolc.h>
46c4762a1bSJed Brown 
47c4762a1bSJed Brown /* (Passive) field for the two variables */
48c4762a1bSJed Brown typedef struct {
49c4762a1bSJed Brown   PetscScalar u,v;
50c4762a1bSJed Brown } Field;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown /* Active field for the two variables */
53c4762a1bSJed Brown typedef struct {
54c4762a1bSJed Brown   adouble u,v;
55c4762a1bSJed Brown } AField;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown /* Application context */
58c4762a1bSJed Brown typedef struct {
59c4762a1bSJed Brown   PetscReal D1,D2,gamma,kappa;
60c4762a1bSJed Brown   AField    **u_a,**f_a;
61c4762a1bSJed Brown   PetscBool aijpc;
62c4762a1bSJed Brown   AdolcCtx  *adctx; /* Automatic differentation support */
63c4762a1bSJed Brown } AppCtx;
64c4762a1bSJed Brown 
65c4762a1bSJed Brown extern PetscErrorCode InitialConditions(DM da,Vec U);
66c4762a1bSJed Brown extern PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y);
67c4762a1bSJed Brown extern PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info,PetscReal t,Field**u,Field**udot,Field**f,void *ptr);
68c4762a1bSJed Brown extern PetscErrorCode IFunctionActive(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr);
69c4762a1bSJed Brown extern PetscErrorCode RHSFunctionActive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr);
70c4762a1bSJed Brown extern PetscErrorCode RHSFunctionPassive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr);
71c4762a1bSJed Brown extern PetscErrorCode IJacobianByHand(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx);
72c4762a1bSJed Brown extern PetscErrorCode IJacobianAdolc(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx);
73c4762a1bSJed Brown extern PetscErrorCode RHSJacobianByHand(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx);
74c4762a1bSJed Brown extern PetscErrorCode RHSJacobianAdolc(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx);
75c4762a1bSJed Brown 
76c4762a1bSJed Brown int main(int argc,char **argv)
77c4762a1bSJed Brown {
78410585f6SHong Zhang   TS             ts;
79410585f6SHong Zhang   Vec            x,r,xdot;
80c4762a1bSJed Brown   DM             da;
81c4762a1bSJed Brown   AppCtx         appctx;
82c4762a1bSJed Brown   AdolcCtx       *adctx;
83c4762a1bSJed Brown   Vec            lambda[1];
84c4762a1bSJed Brown   PetscBool      forwardonly=PETSC_FALSE,implicitform=PETSC_FALSE,byhand=PETSC_FALSE;
85c4762a1bSJed Brown   PetscInt       gxm,gym,i,dofs = 2,ctrl[3] = {0,0,0};
86c4762a1bSJed Brown   PetscScalar    **Seed = NULL,**Rec = NULL,*u_vec;
87c4762a1bSJed Brown   unsigned int   **JP = NULL;
88c4762a1bSJed Brown   ISColoring     iscoloring;
89c4762a1bSJed Brown 
90*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&adctx));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-forwardonly",&forwardonly,NULL));
935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL));
94c4762a1bSJed 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;
955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-aijpc",&appctx.aijpc,NULL));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-no_annotation",&adctx->no_an,NULL));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-jacobian_by_hand",&byhand,NULL));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-adolc_sparse",&adctx->sparse,NULL));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-adolc_sparse_view",&adctx->sparse_view,NULL));
100c4762a1bSJed Brown   appctx.D1    = 8.0e-5;
101c4762a1bSJed Brown   appctx.D2    = 4.0e-5;
102c4762a1bSJed Brown   appctx.gamma = .024;
103c4762a1bSJed Brown   appctx.kappa = .06;
104c4762a1bSJed Brown   appctx.adctx = adctx;
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
108c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,DM_BOUNDARY_PERIODIC,DMDA_STENCIL_STAR,65,65,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,0,"u"));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,1,"v"));
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
117c4762a1bSJed Brown      vectors that are the same types
118c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&x));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&xdot));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124c4762a1bSJed Brown      Create timestepping solver context
125c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSCN));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,da));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
130c4762a1bSJed Brown   if (!implicitform) {
1315f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&appctx));
132c4762a1bSJed Brown   } else {
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDATSSetIFunctionLocal(da,INSERT_VALUES,(DMDATSIFunctionLocal)IFunctionLocalPassive,&appctx));
134c4762a1bSJed Brown   }
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   if (!adctx->no_an) {
1375f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetGhostCorners(da,NULL,NULL,NULL,&gxm,&gym,NULL));
138c4762a1bSJed Brown     adctx->m = dofs*gxm*gym;adctx->n = dofs*gxm*gym; /* Number of dependent and independent variables */
139c4762a1bSJed Brown 
140c4762a1bSJed Brown     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141c4762a1bSJed Brown        Trace function(s) just once
142c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143c4762a1bSJed Brown     if (!implicitform) {
1445f80ce2aSJacob Faibussowitsch       CHKERRQ(RHSFunctionActive(ts,1.0,x,r,&appctx));
145c4762a1bSJed Brown     } else {
1465f80ce2aSJacob Faibussowitsch       CHKERRQ(IFunctionActive(ts,1.0,x,xdot,r,&appctx));
147c4762a1bSJed Brown     }
148c4762a1bSJed Brown 
149c4762a1bSJed Brown     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150c4762a1bSJed Brown       In the case where ADOL-C generates the Jacobian in compressed format,
151c4762a1bSJed Brown       seed and recovery matrices are required. Since the sparsity structure
152c4762a1bSJed Brown       of the Jacobian does not change over the course of the time
153c4762a1bSJed Brown       integration, we can save computational effort by only generating
154c4762a1bSJed Brown       these objects once.
155c4762a1bSJed Brown        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156c4762a1bSJed Brown     if (adctx->sparse) {
157c4762a1bSJed Brown       /*
158c4762a1bSJed Brown          Generate sparsity pattern
159c4762a1bSJed Brown 
160c4762a1bSJed Brown          One way to do this is to use the Jacobian pattern driver `jac_pat`
161c4762a1bSJed Brown          provided by ColPack. Otherwise, it is the responsibility of the user
162c4762a1bSJed Brown          to obtain the sparsity pattern.
163c4762a1bSJed Brown       */
1645f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(adctx->n,&u_vec));
165c4762a1bSJed Brown       JP = (unsigned int **) malloc(adctx->m*sizeof(unsigned int*));
166c4762a1bSJed Brown       jac_pat(1,adctx->m,adctx->n,u_vec,JP,ctrl);
167c4762a1bSJed Brown       if (adctx->sparse_view) {
1685f80ce2aSJacob Faibussowitsch         CHKERRQ(PrintSparsity(MPI_COMM_WORLD,adctx->m,JP));
169c4762a1bSJed Brown       }
170c4762a1bSJed Brown 
171c4762a1bSJed Brown       /*
172c4762a1bSJed Brown         Extract a column colouring
173c4762a1bSJed Brown 
174c4762a1bSJed Brown         For problems using DMDA, colourings can be extracted directly, as
175c4762a1bSJed Brown         follows. There are also ColPack drivers available. Otherwise, it is the
176c4762a1bSJed Brown         responsibility of the user to obtain a suitable colouring.
177c4762a1bSJed Brown       */
1785f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCreateColoring(da,IS_COLORING_LOCAL,&iscoloring));
1795f80ce2aSJacob Faibussowitsch       CHKERRQ(ISColoringGetIS(iscoloring,PETSC_USE_POINTER,&adctx->p,NULL));
180c4762a1bSJed Brown 
181c4762a1bSJed Brown       /* Generate seed matrix to propagate through the forward mode of AD */
1825f80ce2aSJacob Faibussowitsch       CHKERRQ(AdolcMalloc2(adctx->n,adctx->p,&Seed));
1835f80ce2aSJacob Faibussowitsch       CHKERRQ(GenerateSeedMatrix(iscoloring,Seed));
1845f80ce2aSJacob Faibussowitsch       CHKERRQ(ISColoringDestroy(&iscoloring));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown       /*
187c4762a1bSJed Brown         Generate recovery matrix, which is used to recover the Jacobian from
188c4762a1bSJed Brown         compressed format */
1895f80ce2aSJacob Faibussowitsch       CHKERRQ(AdolcMalloc2(adctx->m,adctx->p,&Rec));
1905f80ce2aSJacob Faibussowitsch       CHKERRQ(GetRecoveryMatrix(Seed,JP,adctx->m,adctx->p,Rec));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown       /* Store results and free workspace */
193c4762a1bSJed Brown       adctx->Rec = Rec;
194c4762a1bSJed Brown       for (i=0;i<adctx->m;i++)
195c4762a1bSJed Brown         free(JP[i]);
196c4762a1bSJed Brown       free(JP);
1975f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(u_vec));
198c4762a1bSJed Brown 
199c4762a1bSJed Brown     } else {
200c4762a1bSJed Brown 
201c4762a1bSJed Brown       /*
202c4762a1bSJed Brown         In 'full' Jacobian mode, propagate an identity matrix through the
203c4762a1bSJed Brown         forward mode of AD.
204c4762a1bSJed Brown       */
205c4762a1bSJed Brown       adctx->p = adctx->n;
2065f80ce2aSJacob Faibussowitsch       CHKERRQ(AdolcMalloc2(adctx->n,adctx->p,&Seed));
2075f80ce2aSJacob Faibussowitsch       CHKERRQ(Identity(adctx->n,Seed));
208c4762a1bSJed Brown     }
209c4762a1bSJed Brown     adctx->Seed = Seed;
210c4762a1bSJed Brown   }
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
213c4762a1bSJed Brown      Set Jacobian
214c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
215c4762a1bSJed Brown   if (!implicitform) {
216c4762a1bSJed Brown     if (!byhand) {
2175f80ce2aSJacob Faibussowitsch       CHKERRQ(TSSetRHSJacobian(ts,NULL,NULL,RHSJacobianAdolc,&appctx));
218c4762a1bSJed Brown     } else {
2195f80ce2aSJacob Faibussowitsch       CHKERRQ(TSSetRHSJacobian(ts,NULL,NULL,RHSJacobianByHand,&appctx));
220c4762a1bSJed Brown     }
221c4762a1bSJed Brown   } else {
222c4762a1bSJed Brown     if (appctx.aijpc) {
223c4762a1bSJed Brown       Mat A,B;
224c4762a1bSJed Brown 
2255f80ce2aSJacob Faibussowitsch       CHKERRQ(DMSetMatType(da,MATSELL));
2265f80ce2aSJacob Faibussowitsch       CHKERRQ(DMCreateMatrix(da,&A));
2275f80ce2aSJacob Faibussowitsch       CHKERRQ(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B));
228c4762a1bSJed Brown       /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */
229c4762a1bSJed Brown       if (!byhand) {
2305f80ce2aSJacob Faibussowitsch         CHKERRQ(TSSetIJacobian(ts,A,B,IJacobianAdolc,&appctx));
231c4762a1bSJed Brown       } else {
2325f80ce2aSJacob Faibussowitsch         CHKERRQ(TSSetIJacobian(ts,A,B,IJacobianByHand,&appctx));
233c4762a1bSJed Brown       }
2345f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&A));
2355f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&B));
236c4762a1bSJed Brown     } else {
237c4762a1bSJed Brown       if (!byhand) {
2385f80ce2aSJacob Faibussowitsch         CHKERRQ(TSSetIJacobian(ts,NULL,NULL,IJacobianAdolc,&appctx));
239c4762a1bSJed Brown       } else {
2405f80ce2aSJacob Faibussowitsch         CHKERRQ(TSSetIJacobian(ts,NULL,NULL,IJacobianByHand,&appctx));
241c4762a1bSJed Brown       }
242c4762a1bSJed Brown     }
243c4762a1bSJed Brown   }
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
246c4762a1bSJed Brown      Set initial conditions
247c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(InitialConditions(da,x));
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,x));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252c4762a1bSJed Brown     Have the TS save its trajectory so that TSAdjointSolve() may be used
253c4762a1bSJed Brown     and set solver options
254c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
255c4762a1bSJed Brown   if (!forwardonly) {
2565f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetSaveTrajectory(ts));
2575f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetMaxTime(ts,200.0));
2585f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetTimeStep(ts,0.5));
259c4762a1bSJed Brown   } else {
2605f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetMaxTime(ts,2000.0));
2615f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetTimeStep(ts,10));
262c4762a1bSJed Brown   }
2635f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267c4762a1bSJed Brown      Solve ODE system
268c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2695f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,x));
270c4762a1bSJed Brown   if (!forwardonly) {
271c4762a1bSJed Brown     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272c4762a1bSJed Brown        Start the Adjoint model
273c4762a1bSJed Brown        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2745f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(x,&lambda[0]));
275c4762a1bSJed Brown     /*   Reset initial conditions for the adjoint integration */
2765f80ce2aSJacob Faibussowitsch     CHKERRQ(InitializeLambda(da,lambda[0],0.5,0.5));
2775f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetCostGradients(ts,1,lambda,NULL));
2785f80ce2aSJacob Faibussowitsch     CHKERRQ(TSAdjointSolve(ts));
2795f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&lambda[0]));
280c4762a1bSJed Brown   }
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
283c4762a1bSJed Brown      Free work space.
284c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&xdot));
2865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
2875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
2885f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
289c4762a1bSJed Brown   if (!adctx->no_an) {
2905f80ce2aSJacob Faibussowitsch     if (adctx->sparse) CHKERRQ(AdolcFree2(Rec));
2915f80ce2aSJacob Faibussowitsch     CHKERRQ(AdolcFree2(Seed));
292c4762a1bSJed Brown   }
2935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
2945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(adctx));
295*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
296*b122ec5aSJacob Faibussowitsch   return 0;
297c4762a1bSJed Brown }
298c4762a1bSJed Brown 
299c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U)
300c4762a1bSJed Brown {
301c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
302c4762a1bSJed Brown   Field          **u;
303c4762a1bSJed Brown   PetscReal      hx,hy,x,y;
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   PetscFunctionBegin;
3065f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   hx = 2.5/(PetscReal)Mx;
309c4762a1bSJed Brown   hy = 2.5/(PetscReal)My;
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   /*
312c4762a1bSJed Brown      Get pointers to vector data
313c4762a1bSJed Brown   */
3145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,U,&u));
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   /*
317c4762a1bSJed Brown      Get local grid boundaries
318c4762a1bSJed Brown   */
3195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
320c4762a1bSJed Brown 
321c4762a1bSJed Brown   /*
322c4762a1bSJed Brown      Compute function over the locally owned part of the grid
323c4762a1bSJed Brown   */
324c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
325c4762a1bSJed Brown     y = j*hy;
326c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
327c4762a1bSJed Brown       x = i*hx;
32866baab88SBarry 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;
329c4762a1bSJed Brown       else u[j][i].v = 0.0;
330c4762a1bSJed Brown 
331c4762a1bSJed Brown       u[j][i].u = 1.0 - 2.0*u[j][i].v;
332c4762a1bSJed Brown     }
333c4762a1bSJed Brown   }
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   /*
336c4762a1bSJed Brown      Restore vectors
337c4762a1bSJed Brown   */
3385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,U,&u));
339c4762a1bSJed Brown   PetscFunctionReturn(0);
340c4762a1bSJed Brown }
341c4762a1bSJed Brown 
342c4762a1bSJed Brown PetscErrorCode InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y)
343c4762a1bSJed Brown {
344c4762a1bSJed Brown    PetscInt i,j,Mx,My,xs,ys,xm,ym;
345c4762a1bSJed Brown    Field **l;
346c4762a1bSJed Brown    PetscFunctionBegin;
347c4762a1bSJed Brown 
3485f80ce2aSJacob Faibussowitsch    CHKERRQ(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));
349c4762a1bSJed Brown    /* locate the global i index for x and j index for y */
350c4762a1bSJed Brown    i = (PetscInt)(x*(Mx-1));
351c4762a1bSJed Brown    j = (PetscInt)(y*(My-1));
3525f80ce2aSJacob Faibussowitsch    CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
353c4762a1bSJed Brown 
354c4762a1bSJed Brown    if (xs <= i && i < xs+xm && ys <= j && j < ys+ym) {
355c4762a1bSJed Brown      /* the i,j vertex is on this process */
3565f80ce2aSJacob Faibussowitsch      CHKERRQ(DMDAVecGetArray(da,lambda,&l));
357c4762a1bSJed Brown      l[j][i].u = 1.0;
358c4762a1bSJed Brown      l[j][i].v = 1.0;
3595f80ce2aSJacob Faibussowitsch      CHKERRQ(DMDAVecRestoreArray(da,lambda,&l));
360c4762a1bSJed Brown    }
361c4762a1bSJed Brown    PetscFunctionReturn(0);
362c4762a1bSJed Brown }
363c4762a1bSJed Brown 
364c4762a1bSJed Brown PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info,PetscReal t,Field**u,Field**udot,Field**f,void *ptr)
365c4762a1bSJed Brown {
366c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
367c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym;
368c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
369c4762a1bSJed Brown   PetscScalar    uc,uxx,uyy,vc,vxx,vyy;
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   PetscFunctionBegin;
372c4762a1bSJed Brown   hx = 2.50/(PetscReal)(info->mx); sx = 1.0/(hx*hx);
373c4762a1bSJed Brown   hy = 2.50/(PetscReal)(info->my); sy = 1.0/(hy*hy);
374c4762a1bSJed Brown 
375c4762a1bSJed Brown   /* Get local grid boundaries */
376c4762a1bSJed Brown   xs = info->xs; xm = info->xm; ys = info->ys; ym = info->ym;
377c4762a1bSJed Brown 
378c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
379c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
380c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
381c4762a1bSJed Brown       uc        = u[j][i].u;
382c4762a1bSJed Brown       uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
383c4762a1bSJed Brown       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
384c4762a1bSJed Brown       vc        = u[j][i].v;
385c4762a1bSJed Brown       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
386c4762a1bSJed Brown       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
387c4762a1bSJed Brown       f[j][i].u = udot[j][i].u - appctx->D1*(uxx + uyy) + uc*vc*vc - appctx->gamma*(1.0 - uc);
388c4762a1bSJed Brown       f[j][i].v = udot[j][i].v - appctx->D2*(vxx + vyy) - uc*vc*vc + (appctx->gamma + appctx->kappa)*vc;
389c4762a1bSJed Brown     }
390c4762a1bSJed Brown   }
3915f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(16.0*xm*ym));
392c4762a1bSJed Brown   PetscFunctionReturn(0);
393c4762a1bSJed Brown }
394c4762a1bSJed Brown 
395c4762a1bSJed Brown PetscErrorCode IFunctionActive(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
396c4762a1bSJed Brown {
397c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
398c4762a1bSJed Brown   DM             da;
399c4762a1bSJed Brown   DMDALocalInfo  info;
400c4762a1bSJed Brown   Field          **u,**f,**udot;
401c4762a1bSJed Brown   Vec            localU;
402c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,gxs,gys,gxm,gym;
403c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
404c4762a1bSJed Brown   adouble        uc,uxx,uyy,vc,vxx,vyy;
405c4762a1bSJed Brown   AField         **f_a,*f_c,**u_a,*u_c;
406c4762a1bSJed Brown   PetscScalar    dummy;
407c4762a1bSJed Brown 
408c4762a1bSJed Brown   PetscFunctionBegin;
409c4762a1bSJed Brown 
4105f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
4115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
4125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localU));
413c4762a1bSJed Brown   hx = 2.50/(PetscReal)(info.mx); sx = 1.0/(hx*hx);
414c4762a1bSJed Brown   hy = 2.50/(PetscReal)(info.my); sy = 1.0/(hy*hy);
415c4762a1bSJed Brown   xs = info.xs; xm = info.xm; gxs = info.gxs; gxm = info.gxm;
416c4762a1bSJed Brown   ys = info.ys; ym = info.ym; gys = info.gys; gym = info.gym;
417c4762a1bSJed Brown 
418c4762a1bSJed Brown   /*
419c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
420c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
421c4762a1bSJed Brown      By placing code between these two statements, computations can be
422c4762a1bSJed Brown      done while messages are in transition.
423c4762a1bSJed Brown   */
4245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
4255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
426c4762a1bSJed Brown 
427c4762a1bSJed Brown   /*
428c4762a1bSJed Brown      Get pointers to vector data
429c4762a1bSJed Brown   */
4305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(da,localU,&u));
4315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,F,&f));
4325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(da,Udot,&udot));
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   /*
435c4762a1bSJed Brown     Create contiguous 1-arrays of AFields
436c4762a1bSJed Brown 
437c4762a1bSJed Brown     NOTE: Memory for ADOL-C active variables (such as adouble and AField)
438c4762a1bSJed Brown           cannot be allocated using PetscMalloc, as this does not call the
439c4762a1bSJed Brown           relevant class constructor. Instead, we use the C++ keyword `new`.
440c4762a1bSJed Brown   */
441c4762a1bSJed Brown   u_c = new AField[info.gxm*info.gym];
442c4762a1bSJed Brown   f_c = new AField[info.gxm*info.gym];
443c4762a1bSJed Brown 
444c4762a1bSJed Brown   /* Create corresponding 2-arrays of AFields */
445c4762a1bSJed Brown   u_a = new AField*[info.gym];
446c4762a1bSJed Brown   f_a = new AField*[info.gym];
447c4762a1bSJed Brown 
448c4762a1bSJed Brown   /* Align indices between array types to endow 2d array with ghost points */
4495f80ce2aSJacob Faibussowitsch   CHKERRQ(GiveGhostPoints(da,u_c,&u_a));
4505f80ce2aSJacob Faibussowitsch   CHKERRQ(GiveGhostPoints(da,f_c,&f_a));
451c4762a1bSJed Brown 
452c4762a1bSJed Brown   trace_on(1);  /* Start of active section on tape 1 */
453c4762a1bSJed Brown 
454c4762a1bSJed Brown   /*
455c4762a1bSJed Brown     Mark independence
456c4762a1bSJed Brown 
457c4762a1bSJed Brown     NOTE: Ghost points are marked as independent, in place of the points they represent on
458c4762a1bSJed Brown           other processors / on other boundaries.
459c4762a1bSJed Brown   */
460c4762a1bSJed Brown   for (j=gys; j<gys+gym; j++) {
461c4762a1bSJed Brown     for (i=gxs; i<gxs+gxm; i++) {
462c4762a1bSJed Brown       u_a[j][i].u <<= u[j][i].u;
463c4762a1bSJed Brown       u_a[j][i].v <<= u[j][i].v;
464c4762a1bSJed Brown     }
465c4762a1bSJed Brown   }
466c4762a1bSJed Brown 
467c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
468c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
469c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
470c4762a1bSJed Brown       uc        = u_a[j][i].u;
471c4762a1bSJed Brown       uxx       = (-2.0*uc + u_a[j][i-1].u + u_a[j][i+1].u)*sx;
472c4762a1bSJed Brown       uyy       = (-2.0*uc + u_a[j-1][i].u + u_a[j+1][i].u)*sy;
473c4762a1bSJed Brown       vc        = u_a[j][i].v;
474c4762a1bSJed Brown       vxx       = (-2.0*vc + u_a[j][i-1].v + u_a[j][i+1].v)*sx;
475c4762a1bSJed Brown       vyy       = (-2.0*vc + u_a[j-1][i].v + u_a[j+1][i].v)*sy;
476c4762a1bSJed Brown       f_a[j][i].u = udot[j][i].u - appctx->D1*(uxx + uyy) + uc*vc*vc - appctx->gamma*(1.0 - uc);
477c4762a1bSJed Brown       f_a[j][i].v = udot[j][i].v - appctx->D2*(vxx + vyy) - uc*vc*vc + (appctx->gamma + appctx->kappa)*vc;
478c4762a1bSJed Brown     }
479c4762a1bSJed Brown   }
480c4762a1bSJed Brown 
481c4762a1bSJed Brown   /*
482c4762a1bSJed Brown     Mark dependence
483c4762a1bSJed Brown 
484c4762a1bSJed Brown     NOTE: Marking dependence of dummy variables makes the index notation much simpler when forming
485c4762a1bSJed Brown           the Jacobian later.
486c4762a1bSJed Brown   */
487c4762a1bSJed Brown   for (j=gys; j<gys+gym; j++) {
488c4762a1bSJed Brown     for (i=gxs; i<gxs+gxm; i++) {
489c4762a1bSJed Brown       if ((i < xs) || (i >= xs+xm) || (j < ys) || (j >= ys+ym)) {
490c4762a1bSJed Brown         f_a[j][i].u >>= dummy;
491c4762a1bSJed Brown         f_a[j][i].v >>= dummy;
492c4762a1bSJed Brown       } else {
493c4762a1bSJed Brown         f_a[j][i].u >>= f[j][i].u;
494c4762a1bSJed Brown         f_a[j][i].v >>= f[j][i].v;
495c4762a1bSJed Brown       }
496c4762a1bSJed Brown     }
497c4762a1bSJed Brown   }
498c4762a1bSJed Brown   trace_off();  /* End of active section */
4995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(16.0*xm*ym));
500c4762a1bSJed Brown 
501c4762a1bSJed Brown   /* Restore vectors */
5025f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,F,&f));
5035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(da,localU,&u));
5045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(da,Udot,&udot));
505c4762a1bSJed Brown 
5065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localU));
507410585f6SHong Zhang 
508c4762a1bSJed Brown   /* Destroy AFields appropriately */
509c4762a1bSJed Brown   f_a += info.gys;
510c4762a1bSJed Brown   u_a += info.gys;
511c4762a1bSJed Brown   delete[] f_a;
512c4762a1bSJed Brown   delete[] u_a;
513c4762a1bSJed Brown   delete[] f_c;
514c4762a1bSJed Brown   delete[] u_c;
515c4762a1bSJed Brown 
516c4762a1bSJed Brown   PetscFunctionReturn(0);
517c4762a1bSJed Brown }
518c4762a1bSJed Brown 
519c4762a1bSJed Brown PetscErrorCode RHSFunctionPassive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr)
520c4762a1bSJed Brown {
521c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
522c4762a1bSJed Brown   DM             da;
523c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
524c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
525c4762a1bSJed Brown   PetscScalar    uc,uxx,uyy,vc,vxx,vyy;
526c4762a1bSJed Brown   Field          **u,**f;
527c4762a1bSJed Brown   Vec            localU,localF;
528c4762a1bSJed Brown 
529c4762a1bSJed Brown   PetscFunctionBegin;
5305f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
5315f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
532c4762a1bSJed Brown   hx = 2.50/(PetscReal)(Mx);sx = 1.0/(hx*hx);
533c4762a1bSJed Brown   hy = 2.50/(PetscReal)(My);sy = 1.0/(hy*hy);
5345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localU));
5355f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localF));
536c4762a1bSJed Brown 
537c4762a1bSJed Brown   /*
538c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
539c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
540c4762a1bSJed Brown      By placing code between these two statements, computations can be
541c4762a1bSJed Brown      done while messages are in transition.
542c4762a1bSJed Brown   */
5435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
5445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
5455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F)); // NOTE (1): See (2) below
5465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,F,INSERT_VALUES,localF));
5475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,F,INSERT_VALUES,localF));
548c4762a1bSJed Brown 
549c4762a1bSJed Brown   /*
550c4762a1bSJed Brown      Get pointers to vector data
551c4762a1bSJed Brown   */
5525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(da,localU,&u));
5535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,localF,&f));
554c4762a1bSJed Brown 
555c4762a1bSJed Brown   /*
556c4762a1bSJed Brown      Get local grid boundaries
557c4762a1bSJed Brown   */
5585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
559c4762a1bSJed Brown 
560c4762a1bSJed Brown   /*
561c4762a1bSJed Brown      Compute function over the locally owned part of the grid
562c4762a1bSJed Brown   */
563c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
564c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
565c4762a1bSJed Brown       uc        = u[j][i].u;
566c4762a1bSJed Brown       uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
567c4762a1bSJed Brown       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
568c4762a1bSJed Brown       vc        = u[j][i].v;
569c4762a1bSJed Brown       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
570c4762a1bSJed Brown       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
571c4762a1bSJed Brown       f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);
572c4762a1bSJed Brown       f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc;
573c4762a1bSJed Brown     }
574c4762a1bSJed Brown   }
575c4762a1bSJed Brown 
576c4762a1bSJed Brown   /*
577c4762a1bSJed Brown      Gather global vector, using the 2-step process
578c4762a1bSJed Brown         DMLocalToGlobalBegin(),DMLocalToGlobalEnd().
579c4762a1bSJed Brown 
580c4762a1bSJed Brown      NOTE (2): We need to use ADD_VALUES if boundaries are not of type DM_BOUNDARY_NONE or
581c4762a1bSJed Brown                DM_BOUNDARY_GHOSTED, meaning we should also zero F before addition (see (1) above).
582c4762a1bSJed Brown   */
5835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(da,localF,ADD_VALUES,F));
5845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(da,localF,ADD_VALUES,F));
585c4762a1bSJed Brown 
586c4762a1bSJed Brown   /*
587c4762a1bSJed Brown      Restore vectors
588c4762a1bSJed Brown   */
5895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,localF,&f));
5905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(da,localU,&u));
5915f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localF));
5925f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localU));
5935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(16.0*xm*ym));
594c4762a1bSJed Brown   PetscFunctionReturn(0);
595c4762a1bSJed Brown }
596c4762a1bSJed Brown 
597c4762a1bSJed Brown PetscErrorCode RHSFunctionActive(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr)
598c4762a1bSJed Brown {
599c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
600c4762a1bSJed Brown   DM             da;
601c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,gxs,gys,gxm,gym,Mx,My;
602c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
603c4762a1bSJed Brown   AField         **f_a,*f_c,**u_a,*u_c;
604c4762a1bSJed Brown   adouble        uc,uxx,uyy,vc,vxx,vyy;
605c4762a1bSJed Brown   Field          **u,**f;
606c4762a1bSJed Brown   Vec            localU,localF;
607c4762a1bSJed Brown 
608c4762a1bSJed Brown   PetscFunctionBegin;
6095f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
6105f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
611c4762a1bSJed Brown   hx = 2.50/(PetscReal)(Mx);sx = 1.0/(hx*hx);
612c4762a1bSJed Brown   hy = 2.50/(PetscReal)(My);sy = 1.0/(hy*hy);
6135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localU));
6145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localF));
615c4762a1bSJed Brown 
616c4762a1bSJed Brown   /*
617c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
618c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
619c4762a1bSJed Brown      By placing code between these two statements, computations can be
620c4762a1bSJed Brown      done while messages are in transition.
621c4762a1bSJed Brown   */
6225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
6235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
6245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecZeroEntries(F)); // NOTE (1): See (2) below
6255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,F,INSERT_VALUES,localF));
6265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,F,INSERT_VALUES,localF));
627c4762a1bSJed Brown 
628c4762a1bSJed Brown   /*
629c4762a1bSJed Brown      Get pointers to vector data
630c4762a1bSJed Brown   */
6315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(da,localU,&u));
6325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,localF,&f));
633c4762a1bSJed Brown 
634c4762a1bSJed Brown   /*
635c4762a1bSJed Brown      Get local and ghosted grid boundaries
636c4762a1bSJed Brown   */
6375f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetGhostCorners(da,&gxs,&gys,NULL,&gxm,&gym,NULL));
6385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
639c4762a1bSJed Brown 
640c4762a1bSJed Brown   /*
641c4762a1bSJed Brown     Create contiguous 1-arrays of AFields
642c4762a1bSJed Brown 
643c4762a1bSJed Brown     NOTE: Memory for ADOL-C active variables (such as adouble and AField)
644c4762a1bSJed Brown           cannot be allocated using PetscMalloc, as this does not call the
645c4762a1bSJed Brown           relevant class constructor. Instead, we use the C++ keyword `new`.
646c4762a1bSJed Brown   */
647c4762a1bSJed Brown   u_c = new AField[gxm*gym];
648c4762a1bSJed Brown   f_c = new AField[gxm*gym];
649c4762a1bSJed Brown 
650c4762a1bSJed Brown   /* Create corresponding 2-arrays of AFields */
651c4762a1bSJed Brown   u_a = new AField*[gym];
652c4762a1bSJed Brown   f_a = new AField*[gym];
653c4762a1bSJed Brown 
654c4762a1bSJed Brown   /* Align indices between array types to endow 2d array with ghost points */
6555f80ce2aSJacob Faibussowitsch   CHKERRQ(GiveGhostPoints(da,u_c,&u_a));
6565f80ce2aSJacob Faibussowitsch   CHKERRQ(GiveGhostPoints(da,f_c,&f_a));
657c4762a1bSJed Brown 
658c4762a1bSJed Brown   /*
659c4762a1bSJed Brown      Compute function over the locally owned part of the grid
660c4762a1bSJed Brown   */
661c4762a1bSJed Brown   trace_on(1);  // ----------------------------------------------- Start of active section
662c4762a1bSJed Brown 
663c4762a1bSJed Brown   /*
664c4762a1bSJed Brown     Mark independence
665c4762a1bSJed Brown 
666c4762a1bSJed Brown     NOTE: Ghost points are marked as independent, in place of the points they represent on
667c4762a1bSJed Brown           other processors / on other boundaries.
668c4762a1bSJed Brown   */
669c4762a1bSJed Brown   for (j=gys; j<gys+gym; j++) {
670c4762a1bSJed Brown     for (i=gxs; i<gxs+gxm; i++) {
671c4762a1bSJed Brown       u_a[j][i].u <<= u[j][i].u;
672c4762a1bSJed Brown       u_a[j][i].v <<= u[j][i].v;
673c4762a1bSJed Brown     }
674c4762a1bSJed Brown   }
675c4762a1bSJed Brown 
676c4762a1bSJed Brown   /*
677c4762a1bSJed Brown     Compute function over the locally owned part of the grid
678c4762a1bSJed Brown   */
679c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
680c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
681c4762a1bSJed Brown       uc          = u_a[j][i].u;
682c4762a1bSJed Brown       uxx         = (-2.0*uc + u_a[j][i-1].u + u_a[j][i+1].u)*sx;
683c4762a1bSJed Brown       uyy         = (-2.0*uc + u_a[j-1][i].u + u_a[j+1][i].u)*sy;
684c4762a1bSJed Brown       vc          = u_a[j][i].v;
685c4762a1bSJed Brown       vxx         = (-2.0*vc + u_a[j][i-1].v + u_a[j][i+1].v)*sx;
686c4762a1bSJed Brown       vyy         = (-2.0*vc + u_a[j-1][i].v + u_a[j+1][i].v)*sy;
687c4762a1bSJed Brown       f_a[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);
688c4762a1bSJed Brown       f_a[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc;
689c4762a1bSJed Brown     }
690c4762a1bSJed Brown   }
691c4762a1bSJed Brown   /*
692c4762a1bSJed Brown     Mark dependence
693c4762a1bSJed Brown 
694c4762a1bSJed Brown     NOTE: Ghost points are marked as dependent in order to vastly simplify index notation
695c4762a1bSJed Brown           during Jacobian assembly.
696c4762a1bSJed Brown   */
697c4762a1bSJed Brown   for (j=gys; j<gys+gym; j++) {
698c4762a1bSJed Brown     for (i=gxs; i<gxs+gxm; i++) {
699c4762a1bSJed Brown       f_a[j][i].u >>= f[j][i].u;
700c4762a1bSJed Brown       f_a[j][i].v >>= f[j][i].v;
701c4762a1bSJed Brown     }
702c4762a1bSJed Brown   }
703c4762a1bSJed Brown   trace_off();  // ----------------------------------------------- End of active section
704c4762a1bSJed Brown 
705c4762a1bSJed Brown   /* Test zeroth order scalar evaluation in ADOL-C gives the same result */
706c4762a1bSJed Brown //  if (appctx->adctx->zos) {
7075f80ce2aSJacob Faibussowitsch //    CHKERRQ(TestZOS2d(da,f,u,appctx)); // FIXME: Why does this give nonzero?
708c4762a1bSJed Brown //  }
709c4762a1bSJed Brown 
710c4762a1bSJed Brown   /*
711c4762a1bSJed Brown      Gather global vector, using the 2-step process
712c4762a1bSJed Brown         DMLocalToGlobalBegin(),DMLocalToGlobalEnd().
713c4762a1bSJed Brown 
714c4762a1bSJed Brown      NOTE (2): We need to use ADD_VALUES if boundaries are not of type DM_BOUNDARY_NONE or
715c4762a1bSJed Brown                DM_BOUNDARY_GHOSTED, meaning we should also zero F before addition (see (1) above).
716c4762a1bSJed Brown   */
7175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalBegin(da,localF,ADD_VALUES,F));
7185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMLocalToGlobalEnd(da,localF,ADD_VALUES,F));
719c4762a1bSJed Brown 
720c4762a1bSJed Brown   /*
721c4762a1bSJed Brown      Restore vectors
722c4762a1bSJed Brown   */
7235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,localF,&f));
7245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(da,localU,&u));
7255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localF));
7265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localU));
7275f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(16.0*xm*ym));
728c4762a1bSJed Brown 
729c4762a1bSJed Brown   /* Destroy AFields appropriately */
730c4762a1bSJed Brown   f_a += gys;
731c4762a1bSJed Brown   u_a += gys;
732c4762a1bSJed Brown   delete[] f_a;
733c4762a1bSJed Brown   delete[] u_a;
734c4762a1bSJed Brown   delete[] f_c;
735c4762a1bSJed Brown   delete[] u_c;
736c4762a1bSJed Brown 
737c4762a1bSJed Brown   PetscFunctionReturn(0);
738c4762a1bSJed Brown }
739c4762a1bSJed Brown 
740c4762a1bSJed Brown PetscErrorCode IJacobianAdolc(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
741c4762a1bSJed Brown {
742c4762a1bSJed Brown   AppCtx            *appctx = (AppCtx*)ctx;
743c4762a1bSJed Brown   DM                da;
744a8c08197SHong Zhang   const PetscScalar *u_vec;
745c4762a1bSJed Brown   Vec               localU;
746c4762a1bSJed Brown 
747c4762a1bSJed Brown   PetscFunctionBegin;
7485f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
7495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localU));
750c4762a1bSJed Brown 
751c4762a1bSJed Brown   /*
752c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
753c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
754c4762a1bSJed Brown      By placing code between these two statements, computations can be
755c4762a1bSJed Brown      done while messages are in transition.
756c4762a1bSJed Brown   */
7575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
7585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
759c4762a1bSJed Brown 
760c4762a1bSJed Brown   /* Get pointers to vector data */
7615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(localU,&u_vec));
762c4762a1bSJed Brown 
763c4762a1bSJed Brown   /*
764c4762a1bSJed Brown     Compute Jacobian
765c4762a1bSJed Brown   */
7665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscAdolcComputeIJacobianLocalIDMass(1,A,u_vec,a,appctx->adctx));
767c4762a1bSJed Brown 
768c4762a1bSJed Brown   /*
769c4762a1bSJed Brown      Restore vectors
770c4762a1bSJed Brown   */
7715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(localU,&u_vec));
7725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localU));
773c4762a1bSJed Brown   PetscFunctionReturn(0);
774c4762a1bSJed Brown }
775c4762a1bSJed Brown 
776c4762a1bSJed Brown PetscErrorCode IJacobianByHand(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
777c4762a1bSJed Brown {
778c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
779c4762a1bSJed Brown   DM             da;
780c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
781c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
782c4762a1bSJed Brown   PetscScalar    uc,vc;
783c4762a1bSJed Brown   Field          **u;
784c4762a1bSJed Brown   Vec            localU;
785c4762a1bSJed Brown   MatStencil     stencil[6],rowstencil;
786c4762a1bSJed Brown   PetscScalar    entries[6];
787c4762a1bSJed Brown 
788c4762a1bSJed Brown   PetscFunctionBegin;
7895f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
7905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localU));
7915f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
792c4762a1bSJed Brown 
793c4762a1bSJed Brown   hx = 2.50/(PetscReal)Mx; sx = 1.0/(hx*hx);
794c4762a1bSJed Brown   hy = 2.50/(PetscReal)My; sy = 1.0/(hy*hy);
795c4762a1bSJed Brown 
796c4762a1bSJed Brown   /*
797c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
798c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
799c4762a1bSJed Brown      By placing code between these two statements, computations can be
800c4762a1bSJed Brown      done while messages are in transition.
801c4762a1bSJed Brown   */
8025f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
8035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
804c4762a1bSJed Brown 
805c4762a1bSJed Brown   /*
806c4762a1bSJed Brown      Get pointers to vector data
807c4762a1bSJed Brown   */
8085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(da,localU,&u));
809c4762a1bSJed Brown 
810c4762a1bSJed Brown   /*
811c4762a1bSJed Brown      Get local grid boundaries
812c4762a1bSJed Brown   */
8135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
814c4762a1bSJed Brown 
815c4762a1bSJed Brown   stencil[0].k = 0;
816c4762a1bSJed Brown   stencil[1].k = 0;
817c4762a1bSJed Brown   stencil[2].k = 0;
818c4762a1bSJed Brown   stencil[3].k = 0;
819c4762a1bSJed Brown   stencil[4].k = 0;
820c4762a1bSJed Brown   stencil[5].k = 0;
821c4762a1bSJed Brown   rowstencil.k = 0;
822c4762a1bSJed Brown   /*
823c4762a1bSJed Brown      Compute function over the locally owned part of the grid
824c4762a1bSJed Brown   */
825c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
826c4762a1bSJed Brown 
827c4762a1bSJed Brown     stencil[0].j = j-1;
828c4762a1bSJed Brown     stencil[1].j = j+1;
829c4762a1bSJed Brown     stencil[2].j = j;
830c4762a1bSJed Brown     stencil[3].j = j;
831c4762a1bSJed Brown     stencil[4].j = j;
832c4762a1bSJed Brown     stencil[5].j = j;
833c4762a1bSJed Brown     rowstencil.k = 0; rowstencil.j = j;
834c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
835c4762a1bSJed Brown       uc = u[j][i].u;
836c4762a1bSJed Brown       vc = u[j][i].v;
837c4762a1bSJed Brown 
838c4762a1bSJed Brown       /*      uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
839c4762a1bSJed Brown       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
840c4762a1bSJed Brown 
841c4762a1bSJed Brown       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
842c4762a1bSJed Brown       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
843c4762a1bSJed Brown        f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
844c4762a1bSJed Brown 
845c4762a1bSJed Brown       stencil[0].i = i; stencil[0].c = 0; entries[0] = -appctx->D1*sy;
846c4762a1bSJed Brown       stencil[1].i = i; stencil[1].c = 0; entries[1] = -appctx->D1*sy;
847c4762a1bSJed Brown       stencil[2].i = i-1; stencil[2].c = 0; entries[2] = -appctx->D1*sx;
848c4762a1bSJed Brown       stencil[3].i = i+1; stencil[3].c = 0; entries[3] = -appctx->D1*sx;
849c4762a1bSJed Brown       stencil[4].i = i; stencil[4].c = 0; entries[4] = 2.0*appctx->D1*(sx + sy) + vc*vc + appctx->gamma + a;
850c4762a1bSJed Brown       stencil[5].i = i; stencil[5].c = 1; entries[5] = 2.0*uc*vc;
851c4762a1bSJed Brown       rowstencil.i = i; rowstencil.c = 0;
852c4762a1bSJed Brown 
8535f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
854c4762a1bSJed Brown       if (appctx->aijpc) {
8555f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesStencil(B,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
856c4762a1bSJed Brown       }
857c4762a1bSJed Brown       stencil[0].c = 1; entries[0] = -appctx->D2*sy;
858c4762a1bSJed Brown       stencil[1].c = 1; entries[1] = -appctx->D2*sy;
859c4762a1bSJed Brown       stencil[2].c = 1; entries[2] = -appctx->D2*sx;
860c4762a1bSJed Brown       stencil[3].c = 1; entries[3] = -appctx->D2*sx;
861c4762a1bSJed Brown       stencil[4].c = 1; entries[4] = 2.0*appctx->D2*(sx + sy) - 2.0*uc*vc + appctx->gamma + appctx->kappa + a;
862c4762a1bSJed Brown       stencil[5].c = 0; entries[5] = -vc*vc;
863c4762a1bSJed Brown 
8645f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
865c4762a1bSJed Brown       if (appctx->aijpc) {
8665f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesStencil(B,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
867c4762a1bSJed Brown       }
868c4762a1bSJed Brown       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
869c4762a1bSJed Brown     }
870c4762a1bSJed Brown   }
871c4762a1bSJed Brown 
872c4762a1bSJed Brown   /*
873c4762a1bSJed Brown      Restore vectors
874c4762a1bSJed Brown   */
8755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(19.0*xm*ym));
8765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(da,localU,&u));
8775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localU));
8785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
8795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
8805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
881c4762a1bSJed Brown   if (appctx->aijpc) {
8825f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
8835f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
8845f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
885c4762a1bSJed Brown   }
886c4762a1bSJed Brown   PetscFunctionReturn(0);
887c4762a1bSJed Brown }
888c4762a1bSJed Brown 
889c4762a1bSJed Brown PetscErrorCode RHSJacobianByHand(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
890c4762a1bSJed Brown {
891c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
892c4762a1bSJed Brown   DM             da;
893c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
894c4762a1bSJed Brown   PetscReal      hx,hy,sx,sy;
895c4762a1bSJed Brown   PetscScalar    uc,vc;
896c4762a1bSJed Brown   Field          **u;
897c4762a1bSJed Brown   Vec            localU;
898c4762a1bSJed Brown   MatStencil     stencil[6],rowstencil;
899c4762a1bSJed Brown   PetscScalar    entries[6];
900c4762a1bSJed Brown 
901c4762a1bSJed Brown   PetscFunctionBegin;
9025f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
9035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localU));
9045f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
905c4762a1bSJed Brown 
906c4762a1bSJed Brown   hx = 2.50/(PetscReal)(Mx); sx = 1.0/(hx*hx);
907c4762a1bSJed Brown   hy = 2.50/(PetscReal)(My); sy = 1.0/(hy*hy);
908c4762a1bSJed Brown 
909c4762a1bSJed Brown   /*
910c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
911c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
912c4762a1bSJed Brown      By placing code between these two statements, computations can be
913c4762a1bSJed Brown      done while messages are in transition.
914c4762a1bSJed Brown   */
9155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
9165f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
917c4762a1bSJed Brown 
918c4762a1bSJed Brown   /*
919c4762a1bSJed Brown      Get pointers to vector data
920c4762a1bSJed Brown   */
9215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(da,localU,&u));
922c4762a1bSJed Brown 
923c4762a1bSJed Brown   /*
924c4762a1bSJed Brown      Get local grid boundaries
925c4762a1bSJed Brown   */
9265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
927c4762a1bSJed Brown 
928c4762a1bSJed Brown   stencil[0].k = 0;
929c4762a1bSJed Brown   stencil[1].k = 0;
930c4762a1bSJed Brown   stencil[2].k = 0;
931c4762a1bSJed Brown   stencil[3].k = 0;
932c4762a1bSJed Brown   stencil[4].k = 0;
933c4762a1bSJed Brown   stencil[5].k = 0;
934c4762a1bSJed Brown   rowstencil.k = 0;
935c4762a1bSJed Brown 
936c4762a1bSJed Brown   /*
937c4762a1bSJed Brown      Compute function over the locally owned part of the grid
938c4762a1bSJed Brown   */
939c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
940c4762a1bSJed Brown 
941c4762a1bSJed Brown     stencil[0].j = j-1;
942c4762a1bSJed Brown     stencil[1].j = j+1;
943c4762a1bSJed Brown     stencil[2].j = j;
944c4762a1bSJed Brown     stencil[3].j = j;
945c4762a1bSJed Brown     stencil[4].j = j;
946c4762a1bSJed Brown     stencil[5].j = j;
947c4762a1bSJed Brown     rowstencil.k = 0; rowstencil.j = j;
948c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
949c4762a1bSJed Brown       uc = u[j][i].u;
950c4762a1bSJed Brown       vc = u[j][i].v;
951c4762a1bSJed Brown 
952c4762a1bSJed Brown       /*      uxx       = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
953c4762a1bSJed Brown       uyy       = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
954c4762a1bSJed Brown 
955c4762a1bSJed Brown       vxx       = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
956c4762a1bSJed Brown       vyy       = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
957c4762a1bSJed Brown        f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
958c4762a1bSJed Brown 
959c4762a1bSJed Brown       stencil[0].i = i; stencil[0].c = 0; entries[0] = appctx->D1*sy;
960c4762a1bSJed Brown       stencil[1].i = i; stencil[1].c = 0; entries[1] = appctx->D1*sy;
961c4762a1bSJed Brown       stencil[2].i = i-1; stencil[2].c = 0; entries[2] = appctx->D1*sx;
962c4762a1bSJed Brown       stencil[3].i = i+1; stencil[3].c = 0; entries[3] = appctx->D1*sx;
963c4762a1bSJed Brown       stencil[4].i = i; stencil[4].c = 0; entries[4] = -2.0*appctx->D1*(sx + sy) - vc*vc - appctx->gamma;
964c4762a1bSJed Brown       stencil[5].i = i; stencil[5].c = 1; entries[5] = -2.0*uc*vc;
965c4762a1bSJed Brown       rowstencil.i = i; rowstencil.c = 0;
966c4762a1bSJed Brown 
9675f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
968c4762a1bSJed Brown 
969c4762a1bSJed Brown       stencil[0].c = 1; entries[0] = appctx->D2*sy;
970c4762a1bSJed Brown       stencil[1].c = 1; entries[1] = appctx->D2*sy;
971c4762a1bSJed Brown       stencil[2].c = 1; entries[2] = appctx->D2*sx;
972c4762a1bSJed Brown       stencil[3].c = 1; entries[3] = appctx->D2*sx;
973c4762a1bSJed Brown       stencil[4].c = 1; entries[4] = -2.0*appctx->D2*(sx + sy) + 2.0*uc*vc - appctx->gamma - appctx->kappa;
974c4762a1bSJed Brown       stencil[5].c = 0; entries[5] = vc*vc;
975c4762a1bSJed Brown       rowstencil.c = 1;
976c4762a1bSJed Brown 
9775f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesStencil(A,1,&rowstencil,6,stencil,entries,INSERT_VALUES));
978c4762a1bSJed Brown       /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
979c4762a1bSJed Brown     }
980c4762a1bSJed Brown   }
981c4762a1bSJed Brown 
982c4762a1bSJed Brown   /*
983c4762a1bSJed Brown      Restore vectors
984c4762a1bSJed Brown   */
9855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(19.0*xm*ym));
9865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(da,localU,&u));
9875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localU));
9885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
9895f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
9905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
991c4762a1bSJed Brown   if (appctx->aijpc) {
9925f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
9935f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
9945f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
995c4762a1bSJed Brown   }
996c4762a1bSJed Brown   PetscFunctionReturn(0);
997c4762a1bSJed Brown }
998c4762a1bSJed Brown 
999c4762a1bSJed Brown PetscErrorCode RHSJacobianAdolc(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
1000c4762a1bSJed Brown {
1001c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ctx;
1002c4762a1bSJed Brown   DM             da;
1003c4762a1bSJed Brown   PetscScalar    *u_vec;
1004c4762a1bSJed Brown   Vec            localU;
1005c4762a1bSJed Brown 
1006c4762a1bSJed Brown   PetscFunctionBegin;
10075f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
10085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localU));
1009c4762a1bSJed Brown 
1010c4762a1bSJed Brown   /*
1011c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
1012c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
1013c4762a1bSJed Brown      By placing code between these two statements, computations can be
1014c4762a1bSJed Brown      done while messages are in transition.
1015c4762a1bSJed Brown   */
10165f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
10175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
1018c4762a1bSJed Brown 
1019c4762a1bSJed Brown   /* Get pointers to vector data */
10205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(localU,&u_vec));
1021c4762a1bSJed Brown 
1022c4762a1bSJed Brown   /*
1023c4762a1bSJed Brown     Compute Jacobian
1024c4762a1bSJed Brown   */
10255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscAdolcComputeRHSJacobianLocal(1,A,u_vec,appctx->adctx));
1026c4762a1bSJed Brown 
1027c4762a1bSJed Brown   /*
1028c4762a1bSJed Brown      Restore vectors
1029c4762a1bSJed Brown   */
10305f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(localU,&u_vec));
10315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localU));
1032c4762a1bSJed Brown   PetscFunctionReturn(0);
1033c4762a1bSJed Brown }
1034c4762a1bSJed Brown 
1035c4762a1bSJed Brown /*TEST
1036c4762a1bSJed Brown 
1037c4762a1bSJed Brown   build:
1038c4762a1bSJed Brown     requires: double !complex adolc colpack
1039c4762a1bSJed Brown 
1040c4762a1bSJed Brown   test:
1041c4762a1bSJed Brown     suffix: 1
1042c4762a1bSJed Brown     nsize: 1
1043c4762a1bSJed Brown     args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian
1044c4762a1bSJed Brown     output_file: output/adr_ex5adj_1.out
1045c4762a1bSJed Brown 
1046c4762a1bSJed Brown   test:
1047c4762a1bSJed Brown     suffix: 2
1048c4762a1bSJed Brown     nsize: 1
1049c4762a1bSJed Brown     args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian -implicitform
1050c4762a1bSJed Brown     output_file: output/adr_ex5adj_2.out
1051c4762a1bSJed Brown 
1052c4762a1bSJed Brown   test:
1053c4762a1bSJed Brown     suffix: 3
1054c4762a1bSJed Brown     nsize: 4
1055c4762a1bSJed Brown     args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor
1056c4762a1bSJed Brown     output_file: output/adr_ex5adj_3.out
1057c4762a1bSJed Brown 
1058c4762a1bSJed Brown   test:
1059c4762a1bSJed Brown     suffix: 4
1060c4762a1bSJed Brown     nsize: 4
1061c4762a1bSJed Brown     args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor -implicitform
1062c4762a1bSJed Brown     output_file: output/adr_ex5adj_4.out
1063c4762a1bSJed Brown 
1064c4762a1bSJed Brown   testset:
1065c4762a1bSJed Brown     suffix: 5
1066c4762a1bSJed Brown     nsize: 4
1067c4762a1bSJed Brown     args: -ts_max_steps 10 -da_grid_x 15 -da_grid_y 15 -ts_monitor -ts_adjoint_monitor -adolc_sparse
1068c4762a1bSJed Brown     output_file: output/adr_ex5adj_5.out
1069c4762a1bSJed Brown 
1070c4762a1bSJed Brown   testset:
1071c4762a1bSJed Brown     suffix: 6
1072c4762a1bSJed Brown     nsize: 4
1073c4762a1bSJed Brown     args: -ts_max_steps 10 -da_grid_x 15 -da_grid_y 15 -ts_monitor -ts_adjoint_monitor -adolc_sparse -implicitform
1074c4762a1bSJed Brown     output_file: output/adr_ex5adj_6.out
1075c4762a1bSJed Brown 
1076c4762a1bSJed Brown TEST*/
1077