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.
247addb90fSBarry Smith -aijpc - Set the matrix used to compute the preconditioner 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);
67*2a8381b2SBarry Smith extern PetscErrorCode IJacobianByHand(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, PetscCtx ctx);
68*2a8381b2SBarry Smith extern PetscErrorCode IJacobianAdolc(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, PetscCtx ctx);
69*2a8381b2SBarry Smith extern PetscErrorCode RHSJacobianByHand(TS ts, PetscReal t, Vec U, Mat A, Mat B, PetscCtx ctx);
70*2a8381b2SBarry Smith extern PetscErrorCode RHSJacobianAdolc(TS ts, PetscReal t, Vec U, Mat A, Mat B, PetscCtx ctx);
71c4762a1bSJed Brown
main(int argc,char ** argv)72d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
73d71ae5a4SJacob Faibussowitsch {
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
86327415f7SBarry 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));
919371c9d4SSatish Balay appctx.aijpc = PETSC_FALSE, adctx->no_an = PETSC_FALSE, adctx->sparse = PETSC_FALSE, adctx->sparse_view = PETSC_FALSE;
929371c9d4SSatish Balay adctx->sparse_view_done = PETSC_FALSE;
939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-aijpc", &appctx.aijpc, NULL));
949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_annotation", &adctx->no_an, NULL));
959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-jacobian_by_hand", &byhand, NULL));
969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-adolc_sparse", &adctx->sparse, NULL));
979566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-adolc_sparse_view", &adctx->sparse_view, NULL));
98c4762a1bSJed Brown appctx.D1 = 8.0e-5;
99c4762a1bSJed Brown appctx.D2 = 4.0e-5;
100c4762a1bSJed Brown appctx.gamma = .024;
101c4762a1bSJed Brown appctx.kappa = .06;
102c4762a1bSJed Brown appctx.adctx = adctx;
103c4762a1bSJed Brown
104c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors
106c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1079566063dSJacob 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));
1089566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
1099566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
1109566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "u"));
1119566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "v"));
112c4762a1bSJed Brown
113c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining
115c4762a1bSJed Brown vectors that are the same types
116c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1179566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x));
1189566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r));
1199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xdot));
120c4762a1bSJed Brown
121c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122c4762a1bSJed Brown Create timestepping solver context
123c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1249566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1259566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN));
1269566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da));
1279566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
128c4762a1bSJed Brown if (!implicitform) {
1299566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionPassive, &appctx));
130c4762a1bSJed Brown } else {
1318434afd1SBarry Smith PetscCall(DMDATSSetIFunctionLocal(da, INSERT_VALUES, (DMDATSIFunctionLocalFn *)IFunctionLocalPassive, &appctx));
132c4762a1bSJed Brown }
133c4762a1bSJed Brown
134c4762a1bSJed Brown if (!adctx->no_an) {
1359566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, NULL, NULL, NULL, &gxm, &gym, NULL));
1369371c9d4SSatish Balay adctx->m = dofs * gxm * gym;
1379371c9d4SSatish Balay adctx->n = dofs * gxm * gym; /* Number of dependent and independent variables */
138c4762a1bSJed Brown
139c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140c4762a1bSJed Brown Trace function(s) just once
141c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142c4762a1bSJed Brown if (!implicitform) {
1439566063dSJacob Faibussowitsch PetscCall(RHSFunctionActive(ts, 1.0, x, r, &appctx));
144c4762a1bSJed Brown } else {
1459566063dSJacob Faibussowitsch PetscCall(IFunctionActive(ts, 1.0, x, xdot, r, &appctx));
146c4762a1bSJed Brown }
147c4762a1bSJed Brown
148c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149c4762a1bSJed Brown In the case where ADOL-C generates the Jacobian in compressed format,
150c4762a1bSJed Brown seed and recovery matrices are required. Since the sparsity structure
151c4762a1bSJed Brown of the Jacobian does not change over the course of the time
152c4762a1bSJed Brown integration, we can save computational effort by only generating
153c4762a1bSJed Brown these objects once.
154c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155c4762a1bSJed Brown if (adctx->sparse) {
156c4762a1bSJed Brown /*
157c4762a1bSJed Brown Generate sparsity pattern
158c4762a1bSJed Brown
159c4762a1bSJed Brown One way to do this is to use the Jacobian pattern driver `jac_pat`
160c4762a1bSJed Brown provided by ColPack. Otherwise, it is the responsibility of the user
161c4762a1bSJed Brown to obtain the sparsity pattern.
162c4762a1bSJed Brown */
1639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(adctx->n, &u_vec));
164c4762a1bSJed Brown JP = (unsigned int **)malloc(adctx->m * sizeof(unsigned int *));
165c4762a1bSJed Brown jac_pat(1, adctx->m, adctx->n, u_vec, JP, ctrl);
16648a46eb9SPierre Jolivet if (adctx->sparse_view) PetscCall(PrintSparsity(MPI_COMM_WORLD, adctx->m, JP));
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;
1919371c9d4SSatish Balay for (i = 0; i < adctx->m; i++) free(JP[i]);
192c4762a1bSJed Brown free(JP);
1939566063dSJacob Faibussowitsch PetscCall(PetscFree(u_vec));
194c4762a1bSJed Brown
195c4762a1bSJed Brown } else {
196c4762a1bSJed Brown /*
197c4762a1bSJed Brown In 'full' Jacobian mode, propagate an identity matrix through the
198c4762a1bSJed Brown forward mode of AD.
199c4762a1bSJed Brown */
200c4762a1bSJed Brown adctx->p = adctx->n;
2019566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(adctx->n, adctx->p, &Seed));
2029566063dSJacob Faibussowitsch PetscCall(Identity(adctx->n, Seed));
203c4762a1bSJed Brown }
204c4762a1bSJed Brown adctx->Seed = Seed;
205c4762a1bSJed Brown }
206c4762a1bSJed Brown
207c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208c4762a1bSJed Brown Set Jacobian
209c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
210c4762a1bSJed Brown if (!implicitform) {
211c4762a1bSJed Brown if (!byhand) {
2129566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobianAdolc, &appctx));
213c4762a1bSJed Brown } else {
2149566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobianByHand, &appctx));
215c4762a1bSJed Brown }
216c4762a1bSJed Brown } else {
217c4762a1bSJed Brown if (appctx.aijpc) {
218c4762a1bSJed Brown Mat A, B;
219c4762a1bSJed Brown
2209566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATSELL));
2219566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &A));
2229566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B));
223c4762a1bSJed Brown /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */
224c4762a1bSJed Brown if (!byhand) {
2259566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, B, IJacobianAdolc, &appctx));
226c4762a1bSJed Brown } else {
2279566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, B, IJacobianByHand, &appctx));
228c4762a1bSJed Brown }
2299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
2309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
231c4762a1bSJed Brown } else {
232c4762a1bSJed Brown if (!byhand) {
2339566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobianAdolc, &appctx));
234c4762a1bSJed Brown } else {
2359566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobianByHand, &appctx));
236c4762a1bSJed Brown }
237c4762a1bSJed Brown }
238c4762a1bSJed Brown }
239c4762a1bSJed Brown
240c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
241c4762a1bSJed Brown Set initial conditions
242c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2439566063dSJacob Faibussowitsch PetscCall(InitialConditions(da, x));
2449566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x));
245c4762a1bSJed Brown
246c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
247c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used
248c4762a1bSJed Brown and set solver options
249c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
250c4762a1bSJed Brown if (!forwardonly) {
2519566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts));
2529566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 200.0));
2539566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.5));
254c4762a1bSJed Brown } else {
2559566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 2000.0));
2569566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 10));
257c4762a1bSJed Brown }
2589566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2599566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
260c4762a1bSJed Brown
261c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262c4762a1bSJed Brown Solve ODE system
263c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2649566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x));
265c4762a1bSJed Brown if (!forwardonly) {
266c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267c4762a1bSJed Brown Start the Adjoint model
268c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &lambda[0]));
270c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */
2719566063dSJacob Faibussowitsch PetscCall(InitializeLambda(da, lambda[0], 0.5, 0.5));
2729566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, lambda, NULL));
2739566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts));
2749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0]));
275c4762a1bSJed Brown }
276c4762a1bSJed Brown
277c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
278c4762a1bSJed Brown Free work space.
279c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xdot));
2819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r));
2829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
2839566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
284c4762a1bSJed Brown if (!adctx->no_an) {
2859566063dSJacob Faibussowitsch if (adctx->sparse) PetscCall(AdolcFree2(Rec));
2869566063dSJacob Faibussowitsch PetscCall(AdolcFree2(Seed));
287c4762a1bSJed Brown }
2889566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
2899566063dSJacob Faibussowitsch PetscCall(PetscFree(adctx));
2909566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
291b122ec5aSJacob Faibussowitsch return 0;
292c4762a1bSJed Brown }
293c4762a1bSJed Brown
InitialConditions(DM da,Vec U)294d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(DM da, Vec U)
295d71ae5a4SJacob Faibussowitsch {
296c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My;
297c4762a1bSJed Brown Field **u;
298c4762a1bSJed Brown PetscReal hx, hy, x, y;
299c4762a1bSJed Brown
300c4762a1bSJed Brown PetscFunctionBegin;
3019566063dSJacob 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));
302c4762a1bSJed Brown
303c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx;
304c4762a1bSJed Brown hy = 2.5 / (PetscReal)My;
305c4762a1bSJed Brown
306c4762a1bSJed Brown /*
307c4762a1bSJed Brown Get pointers to vector data
308c4762a1bSJed Brown */
3099566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u));
310c4762a1bSJed Brown
311c4762a1bSJed Brown /*
312c4762a1bSJed Brown Get local grid boundaries
313c4762a1bSJed Brown */
3149566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
315c4762a1bSJed Brown
316c4762a1bSJed Brown /*
317c4762a1bSJed Brown Compute function over the locally owned part of the grid
318c4762a1bSJed Brown */
319c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
320c4762a1bSJed Brown y = j * hy;
321c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
322c4762a1bSJed Brown x = i * hx;
3239371c9d4SSatish Balay if (PetscApproximateGTE(x, 1.0) && PetscApproximateLTE(x, 1.5) && PetscApproximateGTE(y, 1.0) && PetscApproximateLTE(y, 1.5))
3249371c9d4SSatish Balay u[j][i].v = PetscPowReal(PetscSinReal(4.0 * PETSC_PI * x), 2.0) * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * y), 2.0) / 4.0;
325c4762a1bSJed Brown else u[j][i].v = 0.0;
326c4762a1bSJed Brown
327c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v;
328c4762a1bSJed Brown }
329c4762a1bSJed Brown }
330c4762a1bSJed Brown
331c4762a1bSJed Brown /*
332c4762a1bSJed Brown Restore vectors
333c4762a1bSJed Brown */
3349566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u));
3353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
336c4762a1bSJed Brown }
337c4762a1bSJed Brown
InitializeLambda(DM da,Vec lambda,PetscReal x,PetscReal y)338d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeLambda(DM da, Vec lambda, PetscReal x, PetscReal y)
339d71ae5a4SJacob Faibussowitsch {
340c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym;
341c4762a1bSJed Brown Field **l;
342c4762a1bSJed Brown
3434d86920dSPierre Jolivet PetscFunctionBegin;
3449566063dSJacob 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));
345c4762a1bSJed Brown /* locate the global i index for x and j index for y */
346c4762a1bSJed Brown i = (PetscInt)(x * (Mx - 1));
347c4762a1bSJed Brown j = (PetscInt)(y * (My - 1));
3489566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
349c4762a1bSJed Brown
350c4762a1bSJed Brown if (xs <= i && i < xs + xm && ys <= j && j < ys + ym) {
351c4762a1bSJed Brown /* the i,j vertex is on this process */
3529566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, lambda, &l));
353c4762a1bSJed Brown l[j][i].u = 1.0;
354c4762a1bSJed Brown l[j][i].v = 1.0;
3559566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, lambda, &l));
356c4762a1bSJed Brown }
3573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
358c4762a1bSJed Brown }
359c4762a1bSJed Brown
IFunctionLocalPassive(DMDALocalInfo * info,PetscReal t,Field ** u,Field ** udot,Field ** f,void * ptr)360d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info, PetscReal t, Field **u, Field **udot, Field **f, void *ptr)
361d71ae5a4SJacob Faibussowitsch {
362c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ptr;
363c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym;
364c4762a1bSJed Brown PetscReal hx, hy, sx, sy;
365c4762a1bSJed Brown PetscScalar uc, uxx, uyy, vc, vxx, vyy;
366c4762a1bSJed Brown
367c4762a1bSJed Brown PetscFunctionBegin;
368f4f49eeaSPierre Jolivet hx = 2.50 / (PetscReal)info->mx;
3699371c9d4SSatish Balay sx = 1.0 / (hx * hx);
370f4f49eeaSPierre Jolivet hy = 2.50 / (PetscReal)info->my;
3719371c9d4SSatish Balay sy = 1.0 / (hy * hy);
372c4762a1bSJed Brown
373c4762a1bSJed Brown /* Get local grid boundaries */
3749371c9d4SSatish Balay xs = info->xs;
3759371c9d4SSatish Balay xm = info->xm;
3769371c9d4SSatish Balay ys = info->ys;
3779371c9d4SSatish Balay ym = info->ym;
378c4762a1bSJed Brown
379c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
380c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
381c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
382c4762a1bSJed Brown uc = u[j][i].u;
383c4762a1bSJed Brown uxx = (-2.0 * uc + u[j][i - 1].u + u[j][i + 1].u) * sx;
384c4762a1bSJed Brown uyy = (-2.0 * uc + u[j - 1][i].u + u[j + 1][i].u) * sy;
385c4762a1bSJed Brown vc = u[j][i].v;
386c4762a1bSJed Brown vxx = (-2.0 * vc + u[j][i - 1].v + u[j][i + 1].v) * sx;
387c4762a1bSJed Brown vyy = (-2.0 * vc + u[j - 1][i].v + u[j + 1][i].v) * sy;
388c4762a1bSJed Brown f[j][i].u = udot[j][i].u - appctx->D1 * (uxx + uyy) + uc * vc * vc - appctx->gamma * (1.0 - uc);
389c4762a1bSJed Brown f[j][i].v = udot[j][i].v - appctx->D2 * (vxx + vyy) - uc * vc * vc + (appctx->gamma + appctx->kappa) * vc;
390c4762a1bSJed Brown }
391c4762a1bSJed Brown }
3929566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym));
3933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
394c4762a1bSJed Brown }
395c4762a1bSJed Brown
IFunctionActive(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void * ptr)396d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunctionActive(TS ts, PetscReal ftime, Vec U, Vec Udot, Vec F, void *ptr)
397d71ae5a4SJacob Faibussowitsch {
398c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ptr;
399c4762a1bSJed Brown DM da;
400c4762a1bSJed Brown DMDALocalInfo info;
401c4762a1bSJed Brown Field **u, **f, **udot;
402c4762a1bSJed Brown Vec localU;
403c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, gxs, gys, gxm, gym;
404c4762a1bSJed Brown PetscReal hx, hy, sx, sy;
405c4762a1bSJed Brown adouble uc, uxx, uyy, vc, vxx, vyy;
406c4762a1bSJed Brown AField **f_a, *f_c, **u_a, *u_c;
407c4762a1bSJed Brown PetscScalar dummy;
408c4762a1bSJed Brown
409c4762a1bSJed Brown PetscFunctionBegin;
4109566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
4119566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info));
4129566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
413f4f49eeaSPierre Jolivet hx = 2.50 / (PetscReal)info.mx;
4149371c9d4SSatish Balay sx = 1.0 / (hx * hx);
415f4f49eeaSPierre Jolivet hy = 2.50 / (PetscReal)info.my;
4169371c9d4SSatish Balay sy = 1.0 / (hy * hy);
4179371c9d4SSatish Balay xs = info.xs;
4189371c9d4SSatish Balay xm = info.xm;
4199371c9d4SSatish Balay gxs = info.gxs;
4209371c9d4SSatish Balay gxm = info.gxm;
4219371c9d4SSatish Balay ys = info.ys;
4229371c9d4SSatish Balay ym = info.ym;
4239371c9d4SSatish Balay gys = info.gys;
4249371c9d4SSatish Balay gym = info.gym;
425c4762a1bSJed Brown
426c4762a1bSJed Brown /*
427c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
428c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
429c4762a1bSJed Brown By placing code between these two statements, computations can be
430c4762a1bSJed Brown done while messages are in transition.
431c4762a1bSJed Brown */
4329566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
4339566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
434c4762a1bSJed Brown
435c4762a1bSJed Brown /*
436c4762a1bSJed Brown Get pointers to vector data
437c4762a1bSJed Brown */
4389566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u));
4399566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f));
4409566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Udot, &udot));
441c4762a1bSJed Brown
442c4762a1bSJed Brown /*
443c4762a1bSJed Brown Create contiguous 1-arrays of AFields
444c4762a1bSJed Brown
445c4762a1bSJed Brown NOTE: Memory for ADOL-C active variables (such as adouble and AField)
446c4762a1bSJed Brown cannot be allocated using PetscMalloc, as this does not call the
447c4762a1bSJed Brown relevant class constructor. Instead, we use the C++ keyword `new`.
448c4762a1bSJed Brown */
449c4762a1bSJed Brown u_c = new AField[info.gxm * info.gym];
450c4762a1bSJed Brown f_c = new AField[info.gxm * info.gym];
451c4762a1bSJed Brown
452c4762a1bSJed Brown /* Create corresponding 2-arrays of AFields */
453c4762a1bSJed Brown u_a = new AField *[info.gym];
454c4762a1bSJed Brown f_a = new AField *[info.gym];
455c4762a1bSJed Brown
456c4762a1bSJed Brown /* Align indices between array types to endow 2d array with ghost points */
4579566063dSJacob Faibussowitsch PetscCall(GiveGhostPoints(da, u_c, &u_a));
4589566063dSJacob Faibussowitsch PetscCall(GiveGhostPoints(da, f_c, &f_a));
459c4762a1bSJed Brown
460c4762a1bSJed Brown trace_on(1); /* Start of active section on tape 1 */
461c4762a1bSJed Brown
462c4762a1bSJed Brown /*
463c4762a1bSJed Brown Mark independence
464c4762a1bSJed Brown
465c4762a1bSJed Brown NOTE: Ghost points are marked as independent, in place of the points they represent on
466c4762a1bSJed Brown other processors / on other boundaries.
467c4762a1bSJed Brown */
468c4762a1bSJed Brown for (j = gys; j < gys + gym; j++) {
469c4762a1bSJed Brown for (i = gxs; i < gxs + gxm; i++) {
470c4762a1bSJed Brown u_a[j][i].u <<= u[j][i].u;
471c4762a1bSJed Brown u_a[j][i].v <<= u[j][i].v;
472c4762a1bSJed Brown }
473c4762a1bSJed Brown }
474c4762a1bSJed Brown
475c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
476c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
477c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
478c4762a1bSJed Brown uc = u_a[j][i].u;
479c4762a1bSJed Brown uxx = (-2.0 * uc + u_a[j][i - 1].u + u_a[j][i + 1].u) * sx;
480c4762a1bSJed Brown uyy = (-2.0 * uc + u_a[j - 1][i].u + u_a[j + 1][i].u) * sy;
481c4762a1bSJed Brown vc = u_a[j][i].v;
482c4762a1bSJed Brown vxx = (-2.0 * vc + u_a[j][i - 1].v + u_a[j][i + 1].v) * sx;
483c4762a1bSJed Brown vyy = (-2.0 * vc + u_a[j - 1][i].v + u_a[j + 1][i].v) * sy;
484c4762a1bSJed Brown f_a[j][i].u = udot[j][i].u - appctx->D1 * (uxx + uyy) + uc * vc * vc - appctx->gamma * (1.0 - uc);
485c4762a1bSJed Brown f_a[j][i].v = udot[j][i].v - appctx->D2 * (vxx + vyy) - uc * vc * vc + (appctx->gamma + appctx->kappa) * vc;
486c4762a1bSJed Brown }
487c4762a1bSJed Brown }
488c4762a1bSJed Brown
489c4762a1bSJed Brown /*
490c4762a1bSJed Brown Mark dependence
491c4762a1bSJed Brown
492c4762a1bSJed Brown NOTE: Marking dependence of dummy variables makes the index notation much simpler when forming
493c4762a1bSJed Brown the Jacobian later.
494c4762a1bSJed Brown */
495c4762a1bSJed Brown for (j = gys; j < gys + gym; j++) {
496c4762a1bSJed Brown for (i = gxs; i < gxs + gxm; i++) {
497c4762a1bSJed Brown if ((i < xs) || (i >= xs + xm) || (j < ys) || (j >= ys + ym)) {
498c4762a1bSJed Brown f_a[j][i].u >>= dummy;
499c4762a1bSJed Brown f_a[j][i].v >>= dummy;
500c4762a1bSJed Brown } else {
501c4762a1bSJed Brown f_a[j][i].u >>= f[j][i].u;
502c4762a1bSJed Brown f_a[j][i].v >>= f[j][i].v;
503c4762a1bSJed Brown }
504c4762a1bSJed Brown }
505c4762a1bSJed Brown }
506c4762a1bSJed Brown trace_off(); /* End of active section */
5079566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym));
508c4762a1bSJed Brown
509c4762a1bSJed Brown /* Restore vectors */
5109566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f));
5119566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
5129566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Udot, &udot));
513c4762a1bSJed Brown
5149566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
515410585f6SHong Zhang
516c4762a1bSJed Brown /* Destroy AFields appropriately */
517c4762a1bSJed Brown f_a += info.gys;
518c4762a1bSJed Brown u_a += info.gys;
519c4762a1bSJed Brown delete[] f_a;
520c4762a1bSJed Brown delete[] u_a;
521c4762a1bSJed Brown delete[] f_c;
522c4762a1bSJed Brown delete[] u_c;
5233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
524c4762a1bSJed Brown }
525c4762a1bSJed Brown
RHSFunctionPassive(TS ts,PetscReal ftime,Vec U,Vec F,void * ptr)526d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunctionPassive(TS ts, PetscReal ftime, Vec U, Vec F, void *ptr)
527d71ae5a4SJacob Faibussowitsch {
528c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ptr;
529c4762a1bSJed Brown DM da;
530c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My;
531c4762a1bSJed Brown PetscReal hx, hy, sx, sy;
532c4762a1bSJed Brown PetscScalar uc, uxx, uyy, vc, vxx, vyy;
533c4762a1bSJed Brown Field **u, **f;
534c4762a1bSJed Brown Vec localU, localF;
535c4762a1bSJed Brown
536c4762a1bSJed Brown PetscFunctionBegin;
5379566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
5389566063dSJacob 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));
53957508eceSPierre Jolivet hx = 2.50 / (PetscReal)Mx;
5409371c9d4SSatish Balay sx = 1.0 / (hx * hx);
54157508eceSPierre Jolivet hy = 2.50 / (PetscReal)My;
5429371c9d4SSatish Balay sy = 1.0 / (hy * hy);
5439566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
5449566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localF));
545c4762a1bSJed Brown
546c4762a1bSJed Brown /*
547c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
548c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
549c4762a1bSJed Brown By placing code between these two statements, computations can be
550c4762a1bSJed Brown done while messages are in transition.
551c4762a1bSJed Brown */
5529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
5539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
5549566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); // NOTE (1): See (2) below
5559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, F, INSERT_VALUES, localF));
5569566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, F, INSERT_VALUES, localF));
557c4762a1bSJed Brown
558c4762a1bSJed Brown /*
559c4762a1bSJed Brown Get pointers to vector data
560c4762a1bSJed Brown */
5619566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u));
5629566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, localF, &f));
563c4762a1bSJed Brown
564c4762a1bSJed Brown /*
565c4762a1bSJed Brown Get local grid boundaries
566c4762a1bSJed Brown */
5679566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
568c4762a1bSJed Brown
569c4762a1bSJed Brown /*
570c4762a1bSJed Brown Compute function over the locally owned part of the grid
571c4762a1bSJed Brown */
572c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
573c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
574c4762a1bSJed Brown uc = u[j][i].u;
575c4762a1bSJed Brown uxx = (-2.0 * uc + u[j][i - 1].u + u[j][i + 1].u) * sx;
576c4762a1bSJed Brown uyy = (-2.0 * uc + u[j - 1][i].u + u[j + 1][i].u) * sy;
577c4762a1bSJed Brown vc = u[j][i].v;
578c4762a1bSJed Brown vxx = (-2.0 * vc + u[j][i - 1].v + u[j][i + 1].v) * sx;
579c4762a1bSJed Brown vyy = (-2.0 * vc + u[j - 1][i].v + u[j + 1][i].v) * sy;
580c4762a1bSJed Brown f[j][i].u = appctx->D1 * (uxx + uyy) - uc * vc * vc + appctx->gamma * (1.0 - uc);
581c4762a1bSJed Brown f[j][i].v = appctx->D2 * (vxx + vyy) + uc * vc * vc - (appctx->gamma + appctx->kappa) * vc;
582c4762a1bSJed Brown }
583c4762a1bSJed Brown }
584c4762a1bSJed Brown
585c4762a1bSJed Brown /*
586c4762a1bSJed Brown Gather global vector, using the 2-step process
587c4762a1bSJed Brown DMLocalToGlobalBegin(),DMLocalToGlobalEnd().
588c4762a1bSJed Brown
589c4762a1bSJed Brown NOTE (2): We need to use ADD_VALUES if boundaries are not of type DM_BOUNDARY_NONE or
590c4762a1bSJed Brown DM_BOUNDARY_GHOSTED, meaning we should also zero F before addition (see (1) above).
591c4762a1bSJed Brown */
5929566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da, localF, ADD_VALUES, F));
5939566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(da, localF, ADD_VALUES, F));
594c4762a1bSJed Brown
595c4762a1bSJed Brown /*
596c4762a1bSJed Brown Restore vectors
597c4762a1bSJed Brown */
5989566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, localF, &f));
5999566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
6009566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localF));
6019566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
6029566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym));
6033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
604c4762a1bSJed Brown }
605c4762a1bSJed Brown
RHSFunctionActive(TS ts,PetscReal ftime,Vec U,Vec F,void * ptr)606d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunctionActive(TS ts, PetscReal ftime, Vec U, Vec F, void *ptr)
607d71ae5a4SJacob Faibussowitsch {
608c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ptr;
609c4762a1bSJed Brown DM da;
610c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, gxs, gys, gxm, gym, Mx, My;
611c4762a1bSJed Brown PetscReal hx, hy, sx, sy;
612c4762a1bSJed Brown AField **f_a, *f_c, **u_a, *u_c;
613c4762a1bSJed Brown adouble uc, uxx, uyy, vc, vxx, vyy;
614c4762a1bSJed Brown Field **u, **f;
615c4762a1bSJed Brown Vec localU, localF;
616c4762a1bSJed Brown
617c4762a1bSJed Brown PetscFunctionBegin;
6189566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
6199566063dSJacob 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));
62057508eceSPierre Jolivet hx = 2.50 / (PetscReal)Mx;
6219371c9d4SSatish Balay sx = 1.0 / (hx * hx);
62257508eceSPierre Jolivet hy = 2.50 / (PetscReal)My;
6239371c9d4SSatish Balay sy = 1.0 / (hy * hy);
6249566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
6259566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localF));
626c4762a1bSJed Brown
627c4762a1bSJed Brown /*
628c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
629c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
630c4762a1bSJed Brown By placing code between these two statements, computations can be
631c4762a1bSJed Brown done while messages are in transition.
632c4762a1bSJed Brown */
6339566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
6349566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
6359566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); // NOTE (1): See (2) below
6369566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, F, INSERT_VALUES, localF));
6379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, F, INSERT_VALUES, localF));
638c4762a1bSJed Brown
639c4762a1bSJed Brown /*
640c4762a1bSJed Brown Get pointers to vector data
641c4762a1bSJed Brown */
6429566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u));
6439566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, localF, &f));
644c4762a1bSJed Brown
645c4762a1bSJed Brown /*
646c4762a1bSJed Brown Get local and ghosted grid boundaries
647c4762a1bSJed Brown */
6489566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, NULL, &gxm, &gym, NULL));
6499566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
650c4762a1bSJed Brown
651c4762a1bSJed Brown /*
652c4762a1bSJed Brown Create contiguous 1-arrays of AFields
653c4762a1bSJed Brown
654c4762a1bSJed Brown NOTE: Memory for ADOL-C active variables (such as adouble and AField)
655c4762a1bSJed Brown cannot be allocated using PetscMalloc, as this does not call the
656c4762a1bSJed Brown relevant class constructor. Instead, we use the C++ keyword `new`.
657c4762a1bSJed Brown */
658c4762a1bSJed Brown u_c = new AField[gxm * gym];
659c4762a1bSJed Brown f_c = new AField[gxm * gym];
660c4762a1bSJed Brown
661c4762a1bSJed Brown /* Create corresponding 2-arrays of AFields */
662c4762a1bSJed Brown u_a = new AField *[gym];
663c4762a1bSJed Brown f_a = new AField *[gym];
664c4762a1bSJed Brown
665c4762a1bSJed Brown /* Align indices between array types to endow 2d array with ghost points */
6669566063dSJacob Faibussowitsch PetscCall(GiveGhostPoints(da, u_c, &u_a));
6679566063dSJacob Faibussowitsch PetscCall(GiveGhostPoints(da, f_c, &f_a));
668c4762a1bSJed Brown
669c4762a1bSJed Brown /*
670c4762a1bSJed Brown Compute function over the locally owned part of the grid
671c4762a1bSJed Brown */
672c4762a1bSJed Brown trace_on(1); // ----------------------------------------------- Start of active section
673c4762a1bSJed Brown
674c4762a1bSJed Brown /*
675c4762a1bSJed Brown Mark independence
676c4762a1bSJed Brown
677c4762a1bSJed Brown NOTE: Ghost points are marked as independent, in place of the points they represent on
678c4762a1bSJed Brown other processors / on other boundaries.
679c4762a1bSJed Brown */
680c4762a1bSJed Brown for (j = gys; j < gys + gym; j++) {
681c4762a1bSJed Brown for (i = gxs; i < gxs + gxm; i++) {
682c4762a1bSJed Brown u_a[j][i].u <<= u[j][i].u;
683c4762a1bSJed Brown u_a[j][i].v <<= u[j][i].v;
684c4762a1bSJed Brown }
685c4762a1bSJed Brown }
686c4762a1bSJed Brown
687c4762a1bSJed Brown /*
688c4762a1bSJed Brown Compute function over the locally owned part of the grid
689c4762a1bSJed Brown */
690c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
691c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
692c4762a1bSJed Brown uc = u_a[j][i].u;
693c4762a1bSJed Brown uxx = (-2.0 * uc + u_a[j][i - 1].u + u_a[j][i + 1].u) * sx;
694c4762a1bSJed Brown uyy = (-2.0 * uc + u_a[j - 1][i].u + u_a[j + 1][i].u) * sy;
695c4762a1bSJed Brown vc = u_a[j][i].v;
696c4762a1bSJed Brown vxx = (-2.0 * vc + u_a[j][i - 1].v + u_a[j][i + 1].v) * sx;
697c4762a1bSJed Brown vyy = (-2.0 * vc + u_a[j - 1][i].v + u_a[j + 1][i].v) * sy;
698c4762a1bSJed Brown f_a[j][i].u = appctx->D1 * (uxx + uyy) - uc * vc * vc + appctx->gamma * (1.0 - uc);
699c4762a1bSJed Brown f_a[j][i].v = appctx->D2 * (vxx + vyy) + uc * vc * vc - (appctx->gamma + appctx->kappa) * vc;
700c4762a1bSJed Brown }
701c4762a1bSJed Brown }
702c4762a1bSJed Brown /*
703c4762a1bSJed Brown Mark dependence
704c4762a1bSJed Brown
705c4762a1bSJed Brown NOTE: Ghost points are marked as dependent in order to vastly simplify index notation
706c4762a1bSJed Brown during Jacobian assembly.
707c4762a1bSJed Brown */
708c4762a1bSJed Brown for (j = gys; j < gys + gym; j++) {
709c4762a1bSJed Brown for (i = gxs; i < gxs + gxm; i++) {
710c4762a1bSJed Brown f_a[j][i].u >>= f[j][i].u;
711c4762a1bSJed Brown f_a[j][i].v >>= f[j][i].v;
712c4762a1bSJed Brown }
713c4762a1bSJed Brown }
714c4762a1bSJed Brown trace_off(); // ----------------------------------------------- End of active section
715c4762a1bSJed Brown
716c4762a1bSJed Brown /* Test zeroth order scalar evaluation in ADOL-C gives the same result */
717c4762a1bSJed Brown // if (appctx->adctx->zos) {
7189566063dSJacob Faibussowitsch // PetscCall(TestZOS2d(da,f,u,appctx)); // FIXME: Why does this give nonzero?
719c4762a1bSJed Brown // }
720c4762a1bSJed Brown
721c4762a1bSJed Brown /*
722c4762a1bSJed Brown Gather global vector, using the 2-step process
723c4762a1bSJed Brown DMLocalToGlobalBegin(),DMLocalToGlobalEnd().
724c4762a1bSJed Brown
725c4762a1bSJed Brown NOTE (2): We need to use ADD_VALUES if boundaries are not of type DM_BOUNDARY_NONE or
726c4762a1bSJed Brown DM_BOUNDARY_GHOSTED, meaning we should also zero F before addition (see (1) above).
727c4762a1bSJed Brown */
7289566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da, localF, ADD_VALUES, F));
7299566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(da, localF, ADD_VALUES, F));
730c4762a1bSJed Brown
731c4762a1bSJed Brown /*
732c4762a1bSJed Brown Restore vectors
733c4762a1bSJed Brown */
7349566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, localF, &f));
7359566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
7369566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localF));
7379566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
7389566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym));
739c4762a1bSJed Brown
740c4762a1bSJed Brown /* Destroy AFields appropriately */
741c4762a1bSJed Brown f_a += gys;
742c4762a1bSJed Brown u_a += gys;
743c4762a1bSJed Brown delete[] f_a;
744c4762a1bSJed Brown delete[] u_a;
745c4762a1bSJed Brown delete[] f_c;
746c4762a1bSJed Brown delete[] u_c;
7473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
748c4762a1bSJed Brown }
749c4762a1bSJed Brown
IJacobianAdolc(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,PetscCtx ctx)750*2a8381b2SBarry Smith PetscErrorCode IJacobianAdolc(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
751d71ae5a4SJacob Faibussowitsch {
752c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx;
753c4762a1bSJed Brown DM da;
754a8c08197SHong Zhang const PetscScalar *u_vec;
755c4762a1bSJed Brown Vec localU;
756c4762a1bSJed Brown
757c4762a1bSJed Brown PetscFunctionBegin;
7589566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
7599566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
760c4762a1bSJed Brown
761c4762a1bSJed Brown /*
762c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
763c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
764c4762a1bSJed Brown By placing code between these two statements, computations can be
765c4762a1bSJed Brown done while messages are in transition.
766c4762a1bSJed Brown */
7679566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
7689566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
769c4762a1bSJed Brown
770c4762a1bSJed Brown /* Get pointers to vector data */
7719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localU, &u_vec));
772c4762a1bSJed Brown
773c4762a1bSJed Brown /*
774c4762a1bSJed Brown Compute Jacobian
775c4762a1bSJed Brown */
7769566063dSJacob Faibussowitsch PetscCall(PetscAdolcComputeIJacobianLocalIDMass(1, A, u_vec, a, appctx->adctx));
777c4762a1bSJed Brown
778c4762a1bSJed Brown /*
779c4762a1bSJed Brown Restore vectors
780c4762a1bSJed Brown */
7819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localU, &u_vec));
7829566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
7833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
784c4762a1bSJed Brown }
785c4762a1bSJed Brown
IJacobianByHand(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,PetscCtx ctx)786*2a8381b2SBarry Smith PetscErrorCode IJacobianByHand(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
787d71ae5a4SJacob Faibussowitsch {
788c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
789c4762a1bSJed Brown DM da;
790c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym;
791c4762a1bSJed Brown PetscReal hx, hy, sx, sy;
792c4762a1bSJed Brown PetscScalar uc, vc;
793c4762a1bSJed Brown Field **u;
794c4762a1bSJed Brown Vec localU;
795c4762a1bSJed Brown MatStencil stencil[6], rowstencil;
796c4762a1bSJed Brown PetscScalar entries[6];
797c4762a1bSJed Brown
798c4762a1bSJed Brown PetscFunctionBegin;
7999566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
8009566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
8019566063dSJacob 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));
802c4762a1bSJed Brown
8039371c9d4SSatish Balay hx = 2.50 / (PetscReal)Mx;
8049371c9d4SSatish Balay sx = 1.0 / (hx * hx);
8059371c9d4SSatish Balay hy = 2.50 / (PetscReal)My;
8069371c9d4SSatish Balay sy = 1.0 / (hy * hy);
807c4762a1bSJed Brown
808c4762a1bSJed Brown /*
809c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
810c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
811c4762a1bSJed Brown By placing code between these two statements, computations can be
812c4762a1bSJed Brown done while messages are in transition.
813c4762a1bSJed Brown */
8149566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
8159566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
816c4762a1bSJed Brown
817c4762a1bSJed Brown /*
818c4762a1bSJed Brown Get pointers to vector data
819c4762a1bSJed Brown */
8209566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u));
821c4762a1bSJed Brown
822c4762a1bSJed Brown /*
823c4762a1bSJed Brown Get local grid boundaries
824c4762a1bSJed Brown */
8259566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
826c4762a1bSJed Brown
827c4762a1bSJed Brown stencil[0].k = 0;
828c4762a1bSJed Brown stencil[1].k = 0;
829c4762a1bSJed Brown stencil[2].k = 0;
830c4762a1bSJed Brown stencil[3].k = 0;
831c4762a1bSJed Brown stencil[4].k = 0;
832c4762a1bSJed Brown stencil[5].k = 0;
833c4762a1bSJed Brown rowstencil.k = 0;
834c4762a1bSJed Brown /*
835c4762a1bSJed Brown Compute function over the locally owned part of the grid
836c4762a1bSJed Brown */
837c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
838c4762a1bSJed Brown stencil[0].j = j - 1;
839c4762a1bSJed Brown stencil[1].j = j + 1;
840c4762a1bSJed Brown stencil[2].j = j;
841c4762a1bSJed Brown stencil[3].j = j;
842c4762a1bSJed Brown stencil[4].j = j;
843c4762a1bSJed Brown stencil[5].j = j;
8449371c9d4SSatish Balay rowstencil.k = 0;
8459371c9d4SSatish Balay rowstencil.j = j;
846c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
847c4762a1bSJed Brown uc = u[j][i].u;
848c4762a1bSJed Brown vc = u[j][i].v;
849c4762a1bSJed Brown
850c4762a1bSJed Brown /* uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
851c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
852c4762a1bSJed Brown
853c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
854c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
855c4762a1bSJed Brown f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
856c4762a1bSJed Brown
8579371c9d4SSatish Balay stencil[0].i = i;
8589371c9d4SSatish Balay stencil[0].c = 0;
8599371c9d4SSatish Balay entries[0] = -appctx->D1 * sy;
8609371c9d4SSatish Balay stencil[1].i = i;
8619371c9d4SSatish Balay stencil[1].c = 0;
8629371c9d4SSatish Balay entries[1] = -appctx->D1 * sy;
8639371c9d4SSatish Balay stencil[2].i = i - 1;
8649371c9d4SSatish Balay stencil[2].c = 0;
8659371c9d4SSatish Balay entries[2] = -appctx->D1 * sx;
8669371c9d4SSatish Balay stencil[3].i = i + 1;
8679371c9d4SSatish Balay stencil[3].c = 0;
8689371c9d4SSatish Balay entries[3] = -appctx->D1 * sx;
8699371c9d4SSatish Balay stencil[4].i = i;
8709371c9d4SSatish Balay stencil[4].c = 0;
8719371c9d4SSatish Balay entries[4] = 2.0 * appctx->D1 * (sx + sy) + vc * vc + appctx->gamma + a;
8729371c9d4SSatish Balay stencil[5].i = i;
8739371c9d4SSatish Balay stencil[5].c = 1;
8749371c9d4SSatish Balay entries[5] = 2.0 * uc * vc;
8759371c9d4SSatish Balay rowstencil.i = i;
8769371c9d4SSatish Balay rowstencil.c = 0;
877c4762a1bSJed Brown
8789566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
87948a46eb9SPierre Jolivet if (appctx->aijpc) PetscCall(MatSetValuesStencil(B, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
8809371c9d4SSatish Balay stencil[0].c = 1;
8819371c9d4SSatish Balay entries[0] = -appctx->D2 * sy;
8829371c9d4SSatish Balay stencil[1].c = 1;
8839371c9d4SSatish Balay entries[1] = -appctx->D2 * sy;
8849371c9d4SSatish Balay stencil[2].c = 1;
8859371c9d4SSatish Balay entries[2] = -appctx->D2 * sx;
8869371c9d4SSatish Balay stencil[3].c = 1;
8879371c9d4SSatish Balay entries[3] = -appctx->D2 * sx;
8889371c9d4SSatish Balay stencil[4].c = 1;
8899371c9d4SSatish Balay entries[4] = 2.0 * appctx->D2 * (sx + sy) - 2.0 * uc * vc + appctx->gamma + appctx->kappa + a;
8909371c9d4SSatish Balay stencil[5].c = 0;
8919371c9d4SSatish Balay entries[5] = -vc * vc;
892c4762a1bSJed Brown
8939566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
89448a46eb9SPierre Jolivet if (appctx->aijpc) PetscCall(MatSetValuesStencil(B, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
895c4762a1bSJed Brown /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
896c4762a1bSJed Brown }
897c4762a1bSJed Brown }
898c4762a1bSJed Brown
899c4762a1bSJed Brown /*
900c4762a1bSJed Brown Restore vectors
901c4762a1bSJed Brown */
9029566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(19.0 * xm * ym));
9039566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
9049566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
9059566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
9069566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
9079566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
908c4762a1bSJed Brown if (appctx->aijpc) {
9099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
9109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
9119566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
912c4762a1bSJed Brown }
9133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
914c4762a1bSJed Brown }
915c4762a1bSJed Brown
RHSJacobianByHand(TS ts,PetscReal t,Vec U,Mat A,Mat B,PetscCtx ctx)916*2a8381b2SBarry Smith PetscErrorCode RHSJacobianByHand(TS ts, PetscReal t, Vec U, Mat A, Mat B, PetscCtx ctx)
917d71ae5a4SJacob Faibussowitsch {
918c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */
919c4762a1bSJed Brown DM da;
920c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym;
921c4762a1bSJed Brown PetscReal hx, hy, sx, sy;
922c4762a1bSJed Brown PetscScalar uc, vc;
923c4762a1bSJed Brown Field **u;
924c4762a1bSJed Brown Vec localU;
925c4762a1bSJed Brown MatStencil stencil[6], rowstencil;
926c4762a1bSJed Brown PetscScalar entries[6];
927c4762a1bSJed Brown
928c4762a1bSJed Brown PetscFunctionBegin;
9299566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
9309566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
9319566063dSJacob 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));
932c4762a1bSJed Brown
93357508eceSPierre Jolivet hx = 2.50 / (PetscReal)Mx;
9349371c9d4SSatish Balay sx = 1.0 / (hx * hx);
93557508eceSPierre Jolivet hy = 2.50 / (PetscReal)My;
9369371c9d4SSatish Balay sy = 1.0 / (hy * hy);
937c4762a1bSJed Brown
938c4762a1bSJed Brown /*
939c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
940c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
941c4762a1bSJed Brown By placing code between these two statements, computations can be
942c4762a1bSJed Brown done while messages are in transition.
943c4762a1bSJed Brown */
9449566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
9459566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
946c4762a1bSJed Brown
947c4762a1bSJed Brown /*
948c4762a1bSJed Brown Get pointers to vector data
949c4762a1bSJed Brown */
9509566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u));
951c4762a1bSJed Brown
952c4762a1bSJed Brown /*
953c4762a1bSJed Brown Get local grid boundaries
954c4762a1bSJed Brown */
9559566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
956c4762a1bSJed Brown
957c4762a1bSJed Brown stencil[0].k = 0;
958c4762a1bSJed Brown stencil[1].k = 0;
959c4762a1bSJed Brown stencil[2].k = 0;
960c4762a1bSJed Brown stencil[3].k = 0;
961c4762a1bSJed Brown stencil[4].k = 0;
962c4762a1bSJed Brown stencil[5].k = 0;
963c4762a1bSJed Brown rowstencil.k = 0;
964c4762a1bSJed Brown
965c4762a1bSJed Brown /*
966c4762a1bSJed Brown Compute function over the locally owned part of the grid
967c4762a1bSJed Brown */
968c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
969c4762a1bSJed Brown stencil[0].j = j - 1;
970c4762a1bSJed Brown stencil[1].j = j + 1;
971c4762a1bSJed Brown stencil[2].j = j;
972c4762a1bSJed Brown stencil[3].j = j;
973c4762a1bSJed Brown stencil[4].j = j;
974c4762a1bSJed Brown stencil[5].j = j;
9759371c9d4SSatish Balay rowstencil.k = 0;
9769371c9d4SSatish Balay rowstencil.j = j;
977c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
978c4762a1bSJed Brown uc = u[j][i].u;
979c4762a1bSJed Brown vc = u[j][i].v;
980c4762a1bSJed Brown
981c4762a1bSJed Brown /* uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx;
982c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy;
983c4762a1bSJed Brown
984c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx;
985c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy;
986c4762a1bSJed Brown f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/
987c4762a1bSJed Brown
9889371c9d4SSatish Balay stencil[0].i = i;
9899371c9d4SSatish Balay stencil[0].c = 0;
9909371c9d4SSatish Balay entries[0] = appctx->D1 * sy;
9919371c9d4SSatish Balay stencil[1].i = i;
9929371c9d4SSatish Balay stencil[1].c = 0;
9939371c9d4SSatish Balay entries[1] = appctx->D1 * sy;
9949371c9d4SSatish Balay stencil[2].i = i - 1;
9959371c9d4SSatish Balay stencil[2].c = 0;
9969371c9d4SSatish Balay entries[2] = appctx->D1 * sx;
9979371c9d4SSatish Balay stencil[3].i = i + 1;
9989371c9d4SSatish Balay stencil[3].c = 0;
9999371c9d4SSatish Balay entries[3] = appctx->D1 * sx;
10009371c9d4SSatish Balay stencil[4].i = i;
10019371c9d4SSatish Balay stencil[4].c = 0;
10029371c9d4SSatish Balay entries[4] = -2.0 * appctx->D1 * (sx + sy) - vc * vc - appctx->gamma;
10039371c9d4SSatish Balay stencil[5].i = i;
10049371c9d4SSatish Balay stencil[5].c = 1;
10059371c9d4SSatish Balay entries[5] = -2.0 * uc * vc;
10069371c9d4SSatish Balay rowstencil.i = i;
10079371c9d4SSatish Balay rowstencil.c = 0;
1008c4762a1bSJed Brown
10099566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
1010c4762a1bSJed Brown
10119371c9d4SSatish Balay stencil[0].c = 1;
10129371c9d4SSatish Balay entries[0] = appctx->D2 * sy;
10139371c9d4SSatish Balay stencil[1].c = 1;
10149371c9d4SSatish Balay entries[1] = appctx->D2 * sy;
10159371c9d4SSatish Balay stencil[2].c = 1;
10169371c9d4SSatish Balay entries[2] = appctx->D2 * sx;
10179371c9d4SSatish Balay stencil[3].c = 1;
10189371c9d4SSatish Balay entries[3] = appctx->D2 * sx;
10199371c9d4SSatish Balay stencil[4].c = 1;
10209371c9d4SSatish Balay entries[4] = -2.0 * appctx->D2 * (sx + sy) + 2.0 * uc * vc - appctx->gamma - appctx->kappa;
10219371c9d4SSatish Balay stencil[5].c = 0;
10229371c9d4SSatish Balay entries[5] = vc * vc;
1023c4762a1bSJed Brown rowstencil.c = 1;
1024c4762a1bSJed Brown
10259566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES));
1026c4762a1bSJed Brown /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */
1027c4762a1bSJed Brown }
1028c4762a1bSJed Brown }
1029c4762a1bSJed Brown
1030c4762a1bSJed Brown /*
1031c4762a1bSJed Brown Restore vectors
1032c4762a1bSJed Brown */
10339566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(19.0 * xm * ym));
10349566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
10359566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
10369566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
10379566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
10389566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
1039c4762a1bSJed Brown if (appctx->aijpc) {
10409566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
10419566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
10429566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE));
1043c4762a1bSJed Brown }
10443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1045c4762a1bSJed Brown }
1046c4762a1bSJed Brown
RHSJacobianAdolc(TS ts,PetscReal t,Vec U,Mat A,Mat B,PetscCtx ctx)1047*2a8381b2SBarry Smith PetscErrorCode RHSJacobianAdolc(TS ts, PetscReal t, Vec U, Mat A, Mat B, PetscCtx ctx)
1048d71ae5a4SJacob Faibussowitsch {
1049c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx;
1050c4762a1bSJed Brown DM da;
1051c4762a1bSJed Brown PetscScalar *u_vec;
1052c4762a1bSJed Brown Vec localU;
1053c4762a1bSJed Brown
1054c4762a1bSJed Brown PetscFunctionBegin;
10559566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
10569566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU));
1057c4762a1bSJed Brown
1058c4762a1bSJed Brown /*
1059c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process
1060c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
1061c4762a1bSJed Brown By placing code between these two statements, computations can be
1062c4762a1bSJed Brown done while messages are in transition.
1063c4762a1bSJed Brown */
10649566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
10659566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
1066c4762a1bSJed Brown
1067c4762a1bSJed Brown /* Get pointers to vector data */
10689566063dSJacob Faibussowitsch PetscCall(VecGetArray(localU, &u_vec));
1069c4762a1bSJed Brown
1070c4762a1bSJed Brown /*
1071c4762a1bSJed Brown Compute Jacobian
1072c4762a1bSJed Brown */
10739566063dSJacob Faibussowitsch PetscCall(PetscAdolcComputeRHSJacobianLocal(1, A, u_vec, appctx->adctx));
1074c4762a1bSJed Brown
1075c4762a1bSJed Brown /*
1076c4762a1bSJed Brown Restore vectors
1077c4762a1bSJed Brown */
10789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localU, &u_vec));
10799566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU));
10803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1081c4762a1bSJed Brown }
1082c4762a1bSJed Brown
1083c4762a1bSJed Brown /*TEST
1084c4762a1bSJed Brown
1085c4762a1bSJed Brown build:
1086c4762a1bSJed Brown requires: double !complex adolc colpack
1087c4762a1bSJed Brown
1088c4762a1bSJed Brown test:
1089c4762a1bSJed Brown suffix: 1
1090c4762a1bSJed Brown nsize: 1
1091c4762a1bSJed Brown args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian
1092c4762a1bSJed Brown output_file: output/adr_ex5adj_1.out
1093c4762a1bSJed Brown
1094c4762a1bSJed Brown test:
1095c4762a1bSJed Brown suffix: 2
1096c4762a1bSJed Brown nsize: 1
1097c4762a1bSJed Brown args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian -implicitform
1098c4762a1bSJed Brown output_file: output/adr_ex5adj_2.out
1099c4762a1bSJed Brown
1100c4762a1bSJed Brown test:
1101c4762a1bSJed Brown suffix: 3
1102c4762a1bSJed Brown nsize: 4
1103c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor
1104c4762a1bSJed Brown output_file: output/adr_ex5adj_3.out
1105c4762a1bSJed Brown
1106c4762a1bSJed Brown test:
1107c4762a1bSJed Brown suffix: 4
1108c4762a1bSJed Brown nsize: 4
1109c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor -implicitform
1110c4762a1bSJed Brown output_file: output/adr_ex5adj_4.out
1111c4762a1bSJed Brown
1112c4762a1bSJed Brown testset:
1113c4762a1bSJed Brown suffix: 5
1114c4762a1bSJed Brown nsize: 4
1115c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 15 -da_grid_y 15 -ts_monitor -ts_adjoint_monitor -adolc_sparse
1116c4762a1bSJed Brown output_file: output/adr_ex5adj_5.out
1117c4762a1bSJed Brown
1118c4762a1bSJed Brown testset:
1119c4762a1bSJed Brown suffix: 6
1120c4762a1bSJed Brown nsize: 4
1121c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 15 -da_grid_y 15 -ts_monitor -ts_adjoint_monitor -adolc_sparse -implicitform
1122c4762a1bSJed Brown output_file: output/adr_ex5adj_6.out
1123c4762a1bSJed Brown
1124c4762a1bSJed Brown TEST*/
1125