1c4762a1bSJed Brown static char help[] = "Demonstrates automatic, matrix-free Jacobian generation using ADOL-C for a time-dependent PDE in 2d, solved using implicit timestepping.\n"; 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 /* ------------------------------------------------------------------------ 10c4762a1bSJed Brown See ../advection-diffusion-reaction/ex5 for a description of the problem 11c4762a1bSJed Brown ------------------------------------------------------------------------- */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown #include <petscdmda.h> 14c4762a1bSJed Brown #include <petscts.h> 15c4762a1bSJed Brown #include "adolc-utils/init.cxx" 16c4762a1bSJed Brown #include "adolc-utils/matfree.cxx" 17c4762a1bSJed Brown #include <adolc/adolc.h> 18c4762a1bSJed Brown 19c4762a1bSJed Brown /* (Passive) field for the two variables */ 20c4762a1bSJed Brown typedef struct { 21c4762a1bSJed Brown PetscScalar u, v; 22c4762a1bSJed Brown } Field; 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* Active field for the two variables */ 25c4762a1bSJed Brown typedef struct { 26c4762a1bSJed Brown adouble u, v; 27c4762a1bSJed Brown } AField; 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* Application context */ 30c4762a1bSJed Brown typedef struct { 31c4762a1bSJed Brown PetscReal D1, D2, gamma, kappa; 32c4762a1bSJed Brown AField **u_a, **f_a; 33c4762a1bSJed Brown AdolcCtx *adctx; /* Automatic differentation support */ 34c4762a1bSJed Brown } AppCtx; 35c4762a1bSJed Brown 36c4762a1bSJed Brown extern PetscErrorCode InitialConditions(DM da, Vec U); 37c4762a1bSJed Brown extern PetscErrorCode InitializeLambda(DM da, Vec lambda, PetscReal x, PetscReal y); 38c4762a1bSJed Brown extern PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info, PetscReal t, Field **u, Field **udot, Field **f, void *ptr); 39c4762a1bSJed Brown extern PetscErrorCode IFunctionActive(TS ts, PetscReal ftime, Vec U, Vec Udot, Vec F, void *ptr); 40c4762a1bSJed Brown extern PetscErrorCode IJacobianMatFree(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A_shell, Mat B, void *ctx); 41c4762a1bSJed Brown 42d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 43d71ae5a4SJacob Faibussowitsch { 44c4762a1bSJed Brown TS ts; /* ODE integrator */ 45c4762a1bSJed Brown Vec x, r; /* solution, residual */ 46c4762a1bSJed Brown DM da; 47c4762a1bSJed Brown AppCtx appctx; /* Application context */ 48c4762a1bSJed Brown AdolcMatCtx matctx; /* Matrix (free) context */ 49c4762a1bSJed Brown Vec lambda[1]; 50c4762a1bSJed Brown PetscBool forwardonly = PETSC_FALSE; 51c4762a1bSJed Brown Mat A; /* (Matrix free) Jacobian matrix */ 52c4762a1bSJed Brown PetscInt gxm, gym; 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 55c4762a1bSJed Brown Initialize program 56c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 57327415f7SBarry Smith PetscFunctionBeginUser; 589566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-forwardonly", &forwardonly, NULL)); 60c4762a1bSJed Brown PetscFunctionBeginUser; 61c4762a1bSJed Brown appctx.D1 = 8.0e-5; 62c4762a1bSJed Brown appctx.D2 = 4.0e-5; 63c4762a1bSJed Brown appctx.gamma = .024; 64c4762a1bSJed Brown appctx.kappa = .06; 659566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("df/dx forward", MAT_CLASSID, &matctx.event1)); 669566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("df/d(xdot) forward", MAT_CLASSID, &matctx.event2)); 679566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("df/dx reverse", MAT_CLASSID, &matctx.event3)); 689566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("df/d(xdot) reverse", MAT_CLASSID, &matctx.event4)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 71c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 72c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 739566063dSJacob 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)); 749566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 759566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 769566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "u")); 779566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "v")); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 80c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 81c4762a1bSJed Brown vectors that are the same types 82c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 839566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 87c4762a1bSJed Brown Create matrix free context and specify usage of PETSc-ADOL-C drivers 88c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 899566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATSHELL)); 909566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &A)); 919566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A, &matctx)); 929566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_MULT, (void (*)(void))PetscAdolcIJacobianVectorProductIDMass)); 939566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A, MATOP_MULT_TRANSPOSE, (void (*)(void))PetscAdolcIJacobianTransposeVectorProductIDMass)); 949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &matctx.X)); 959566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &matctx.Xdot)); 969566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &matctx.localX0)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 99c4762a1bSJed Brown Create timestepping solver context 100c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1019566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1029566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN)); 1039566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 1049566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1059566063dSJacob Faibussowitsch PetscCall(DMDATSSetIFunctionLocal(da, INSERT_VALUES, (DMDATSIFunctionLocal)IFunctionLocalPassive, &appctx)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 108c4762a1bSJed Brown Some data required for matrix-free context 109c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1109566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, NULL, NULL, NULL, &gxm, &gym, NULL)); 1119371c9d4SSatish Balay matctx.m = 2 * gxm * gym; 1129371c9d4SSatish Balay matctx.n = 2 * gxm * gym; /* Number of dependent and independent variables */ 113c4762a1bSJed Brown matctx.flg = PETSC_FALSE; /* Flag for reverse mode */ 114c4762a1bSJed Brown matctx.tag1 = 1; /* Tape identifier */ 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 117c4762a1bSJed Brown Trace function just once 118c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1199566063dSJacob Faibussowitsch PetscCall(PetscNew(&appctx.adctx)); 1209566063dSJacob Faibussowitsch PetscCall(IFunctionActive(ts, 1., x, matctx.Xdot, r, &appctx)); 1219566063dSJacob Faibussowitsch PetscCall(PetscFree(appctx.adctx)); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 124c4762a1bSJed Brown Set Jacobian. In this case, IJacobian simply acts to pass context 125c4762a1bSJed Brown information to the matrix-free Jacobian vector product. 126c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1279566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, IJacobianMatFree, &appctx)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 130c4762a1bSJed Brown Set initial conditions 131c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1329566063dSJacob Faibussowitsch PetscCall(InitialConditions(da, x)); 1339566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used 137c4762a1bSJed Brown and set solver options 138c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 139c4762a1bSJed Brown if (!forwardonly) { 1409566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 1419566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 200.0)); 1429566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.5)); 143c4762a1bSJed Brown } else { 1449566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 2000.0)); 1459566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 10)); 146c4762a1bSJed Brown } 1479566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 1489566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 151c4762a1bSJed Brown Solve ODE system 152c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1539566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 154c4762a1bSJed Brown if (!forwardonly) { 155c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156c4762a1bSJed Brown Start the Adjoint model 157c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &lambda[0])); 159c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 1609566063dSJacob Faibussowitsch PetscCall(InitializeLambda(da, lambda[0], 0.5, 0.5)); 1619566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, lambda, NULL)); 1629566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 1639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0])); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 167c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 168c4762a1bSJed Brown are no longer needed. 169c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1709566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &matctx.localX0)); 1719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 1729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&matctx.X)); 1739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&matctx.Xdot)); 1749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1769566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1779566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 178c4762a1bSJed Brown 1799566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 180b122ec5aSJacob Faibussowitsch return 0; 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(DM da, Vec U) 184d71ae5a4SJacob Faibussowitsch { 185c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My; 186c4762a1bSJed Brown Field **u; 187c4762a1bSJed Brown PetscReal hx, hy, x, y; 188c4762a1bSJed Brown 189c4762a1bSJed Brown PetscFunctionBegin; 1909566063dSJacob 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)); 191c4762a1bSJed Brown 192c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx; 193c4762a1bSJed Brown hy = 2.5 / (PetscReal)My; 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* 196c4762a1bSJed Brown Get pointers to vector data 197c4762a1bSJed Brown */ 1989566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 199c4762a1bSJed Brown 200c4762a1bSJed Brown /* 201c4762a1bSJed Brown Get local grid boundaries 202c4762a1bSJed Brown */ 2039566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* 206c4762a1bSJed Brown Compute function over the locally owned part of the grid 207c4762a1bSJed Brown */ 208c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 209c4762a1bSJed Brown y = j * hy; 210c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 211c4762a1bSJed Brown x = i * hx; 2129371c9d4SSatish Balay if (PetscApproximateGTE(x, 1.0) && PetscApproximateLTE(x, 1.5) && PetscApproximateGTE(y, 1.0) && PetscApproximateLTE(y, 1.5)) 2139371c9d4SSatish 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; 214c4762a1bSJed Brown else u[j][i].v = 0.0; 215c4762a1bSJed Brown 216c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v; 217c4762a1bSJed Brown } 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* 221c4762a1bSJed Brown Restore vectors 222c4762a1bSJed Brown */ 2239566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 224*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 225c4762a1bSJed Brown } 226c4762a1bSJed Brown 227d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeLambda(DM da, Vec lambda, PetscReal x, PetscReal y) 228d71ae5a4SJacob Faibussowitsch { 229c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym; 230c4762a1bSJed Brown Field **l; 231c4762a1bSJed Brown 232410585f6SHong Zhang PetscFunctionBegin; 2339566063dSJacob 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)); 234c4762a1bSJed Brown /* locate the global i index for x and j index for y */ 235c4762a1bSJed Brown i = (PetscInt)(x * (Mx - 1)); 236c4762a1bSJed Brown j = (PetscInt)(y * (My - 1)); 2379566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 238c4762a1bSJed Brown 239c4762a1bSJed Brown if (xs <= i && i < xs + xm && ys <= j && j < ys + ym) { 240c4762a1bSJed Brown /* the i,j vertex is on this process */ 2419566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, lambda, &l)); 242c4762a1bSJed Brown l[j][i].u = 1.0; 243c4762a1bSJed Brown l[j][i].v = 1.0; 2449566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, lambda, &l)); 245c4762a1bSJed Brown } 246*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 247c4762a1bSJed Brown } 248c4762a1bSJed Brown 249d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info, PetscReal t, Field **u, Field **udot, Field **f, void *ptr) 250d71ae5a4SJacob Faibussowitsch { 251c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ptr; 252c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym; 253c4762a1bSJed Brown PetscReal hx, hy, sx, sy; 254c4762a1bSJed Brown PetscScalar uc, uxx, uyy, vc, vxx, vyy; 255c4762a1bSJed Brown 256c4762a1bSJed Brown PetscFunctionBegin; 2579371c9d4SSatish Balay hx = 2.50 / (PetscReal)(info->mx); 2589371c9d4SSatish Balay sx = 1.0 / (hx * hx); 2599371c9d4SSatish Balay hy = 2.50 / (PetscReal)(info->my); 2609371c9d4SSatish Balay sy = 1.0 / (hy * hy); 261c4762a1bSJed Brown 262c4762a1bSJed Brown /* Get local grid boundaries */ 2639371c9d4SSatish Balay xs = info->xs; 2649371c9d4SSatish Balay xm = info->xm; 2659371c9d4SSatish Balay ys = info->ys; 2669371c9d4SSatish Balay ym = info->ym; 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 269c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 270c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 271c4762a1bSJed Brown uc = u[j][i].u; 272c4762a1bSJed Brown uxx = (-2.0 * uc + u[j][i - 1].u + u[j][i + 1].u) * sx; 273c4762a1bSJed Brown uyy = (-2.0 * uc + u[j - 1][i].u + u[j + 1][i].u) * sy; 274c4762a1bSJed Brown vc = u[j][i].v; 275c4762a1bSJed Brown vxx = (-2.0 * vc + u[j][i - 1].v + u[j][i + 1].v) * sx; 276c4762a1bSJed Brown vyy = (-2.0 * vc + u[j - 1][i].v + u[j + 1][i].v) * sy; 277c4762a1bSJed Brown f[j][i].u = udot[j][i].u - appctx->D1 * (uxx + uyy) + uc * vc * vc - appctx->gamma * (1.0 - uc); 278c4762a1bSJed Brown f[j][i].v = udot[j][i].v - appctx->D2 * (vxx + vyy) - uc * vc * vc + (appctx->gamma + appctx->kappa) * vc; 279c4762a1bSJed Brown } 280c4762a1bSJed Brown } 2819566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym)); 282*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 283c4762a1bSJed Brown } 284c4762a1bSJed Brown 285d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunctionActive(TS ts, PetscReal ftime, Vec U, Vec Udot, Vec F, void *ptr) 286d71ae5a4SJacob Faibussowitsch { 287c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ptr; 288c4762a1bSJed Brown DM da; 289c4762a1bSJed Brown DMDALocalInfo info; 290c4762a1bSJed Brown Field **u, **f, **udot; 291c4762a1bSJed Brown Vec localU; 292c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, gxs, gys, gxm, gym; 293c4762a1bSJed Brown PetscReal hx, hy, sx, sy; 294c4762a1bSJed Brown adouble uc, uxx, uyy, vc, vxx, vyy; 295c4762a1bSJed Brown AField **f_a, *f_c, **u_a, *u_c; 296c4762a1bSJed Brown PetscScalar dummy; 297c4762a1bSJed Brown 298c4762a1bSJed Brown PetscFunctionBegin; 2999566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 3009566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 3019566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 3029371c9d4SSatish Balay hx = 2.50 / (PetscReal)(info.mx); 3039371c9d4SSatish Balay sx = 1.0 / (hx * hx); 3049371c9d4SSatish Balay hy = 2.50 / (PetscReal)(info.my); 3059371c9d4SSatish Balay sy = 1.0 / (hy * hy); 3069371c9d4SSatish Balay xs = info.xs; 3079371c9d4SSatish Balay xm = info.xm; 3089371c9d4SSatish Balay gxs = info.gxs; 3099371c9d4SSatish Balay gxm = info.gxm; 3109371c9d4SSatish Balay ys = info.ys; 3119371c9d4SSatish Balay ym = info.ym; 3129371c9d4SSatish Balay gys = info.gys; 3139371c9d4SSatish Balay gym = info.gym; 314c4762a1bSJed Brown 315c4762a1bSJed Brown /* 316c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 317c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 318c4762a1bSJed Brown By placing code between these two statements, computations can be 319c4762a1bSJed Brown done while messages are in transition. 320c4762a1bSJed Brown */ 3219566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 3229566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 323c4762a1bSJed Brown 324c4762a1bSJed Brown /* 325c4762a1bSJed Brown Get pointers to vector data 326c4762a1bSJed Brown */ 3279566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u)); 3289566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 3299566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Udot, &udot)); 330c4762a1bSJed Brown 331c4762a1bSJed Brown /* 332c4762a1bSJed Brown Create contiguous 1-arrays of AFields 333c4762a1bSJed Brown 334c4762a1bSJed Brown NOTE: Memory for ADOL-C active variables (such as adouble and AField) 335c4762a1bSJed Brown cannot be allocated using PetscMalloc, as this does not call the 336c4762a1bSJed Brown relevant class constructor. Instead, we use the C++ keyword `new`. 337c4762a1bSJed Brown */ 338c4762a1bSJed Brown u_c = new AField[info.gxm * info.gym]; 339c4762a1bSJed Brown f_c = new AField[info.gxm * info.gym]; 340c4762a1bSJed Brown 341c4762a1bSJed Brown /* Create corresponding 2-arrays of AFields */ 342c4762a1bSJed Brown u_a = new AField *[info.gym]; 343c4762a1bSJed Brown f_a = new AField *[info.gym]; 344c4762a1bSJed Brown 345c4762a1bSJed Brown /* Align indices between array types to endow 2d array with ghost points */ 3469566063dSJacob Faibussowitsch PetscCall(GiveGhostPoints(da, u_c, &u_a)); 3479566063dSJacob Faibussowitsch PetscCall(GiveGhostPoints(da, f_c, &f_a)); 348c4762a1bSJed Brown 349c4762a1bSJed Brown trace_on(1); /* Start of active section on tape 1 */ 350c4762a1bSJed Brown 351c4762a1bSJed Brown /* 352c4762a1bSJed Brown Mark independence 353c4762a1bSJed Brown 354c4762a1bSJed Brown NOTE: Ghost points are marked as independent, in place of the points they represent on 355c4762a1bSJed Brown other processors / on other boundaries. 356c4762a1bSJed Brown */ 357c4762a1bSJed Brown for (j = gys; j < gys + gym; j++) { 358c4762a1bSJed Brown for (i = gxs; i < gxs + gxm; i++) { 359c4762a1bSJed Brown u_a[j][i].u <<= u[j][i].u; 360c4762a1bSJed Brown u_a[j][i].v <<= u[j][i].v; 361c4762a1bSJed Brown } 362c4762a1bSJed Brown } 363c4762a1bSJed Brown 364c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 365c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 366c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 367c4762a1bSJed Brown uc = u_a[j][i].u; 368c4762a1bSJed Brown uxx = (-2.0 * uc + u_a[j][i - 1].u + u_a[j][i + 1].u) * sx; 369c4762a1bSJed Brown uyy = (-2.0 * uc + u_a[j - 1][i].u + u_a[j + 1][i].u) * sy; 370c4762a1bSJed Brown vc = u_a[j][i].v; 371c4762a1bSJed Brown vxx = (-2.0 * vc + u_a[j][i - 1].v + u_a[j][i + 1].v) * sx; 372c4762a1bSJed Brown vyy = (-2.0 * vc + u_a[j - 1][i].v + u_a[j + 1][i].v) * sy; 373c4762a1bSJed Brown f_a[j][i].u = udot[j][i].u - appctx->D1 * (uxx + uyy) + uc * vc * vc - appctx->gamma * (1.0 - uc); 374c4762a1bSJed Brown f_a[j][i].v = udot[j][i].v - appctx->D2 * (vxx + vyy) - uc * vc * vc + (appctx->gamma + appctx->kappa) * vc; 375c4762a1bSJed Brown } 376c4762a1bSJed Brown } 377c4762a1bSJed Brown 378c4762a1bSJed Brown /* 379c4762a1bSJed Brown Mark dependence 380c4762a1bSJed Brown 381c4762a1bSJed Brown NOTE: Marking dependence of dummy variables makes the index notation much simpler when forming 382c4762a1bSJed Brown the Jacobian later. 383c4762a1bSJed Brown */ 384c4762a1bSJed Brown for (j = gys; j < gys + gym; j++) { 385c4762a1bSJed Brown for (i = gxs; i < gxs + gxm; i++) { 386c4762a1bSJed Brown if ((i < xs) || (i >= xs + xm) || (j < ys) || (j >= ys + ym)) { 387c4762a1bSJed Brown f_a[j][i].u >>= dummy; 388c4762a1bSJed Brown f_a[j][i].v >>= dummy; 389c4762a1bSJed Brown } else { 390c4762a1bSJed Brown f_a[j][i].u >>= f[j][i].u; 391c4762a1bSJed Brown f_a[j][i].v >>= f[j][i].v; 392c4762a1bSJed Brown } 393c4762a1bSJed Brown } 394c4762a1bSJed Brown } 395c4762a1bSJed Brown trace_off(); /* End of active section */ 3969566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym)); 397c4762a1bSJed Brown 398c4762a1bSJed Brown /* Restore vectors */ 3999566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 4009566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u)); 4019566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Udot, &udot)); 402c4762a1bSJed Brown 4039566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 404410585f6SHong Zhang 405c4762a1bSJed Brown /* Destroy AFields appropriately */ 406c4762a1bSJed Brown f_a += info.gys; 407c4762a1bSJed Brown u_a += info.gys; 408c4762a1bSJed Brown delete[] f_a; 409c4762a1bSJed Brown delete[] u_a; 410c4762a1bSJed Brown delete[] f_c; 411c4762a1bSJed Brown delete[] u_c; 412*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 413c4762a1bSJed Brown } 414c4762a1bSJed Brown 415c4762a1bSJed Brown /* 416c4762a1bSJed Brown Simply acts to pass TS information to the AdolcMatCtx 417c4762a1bSJed Brown */ 418d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobianMatFree(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A_shell, Mat B, void *ctx) 419d71ae5a4SJacob Faibussowitsch { 420c4762a1bSJed Brown AdolcMatCtx *mctx; 421c4762a1bSJed Brown DM da; 422c4762a1bSJed Brown 423c4762a1bSJed Brown PetscFunctionBeginUser; 4249566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A_shell, &mctx)); 425c4762a1bSJed Brown 426c4762a1bSJed Brown mctx->time = t; 427c4762a1bSJed Brown mctx->shift = a; 428c4762a1bSJed Brown if (mctx->ts != ts) mctx->ts = ts; 4299566063dSJacob Faibussowitsch PetscCall(VecCopy(X, mctx->X)); 4309566063dSJacob Faibussowitsch PetscCall(VecCopy(Xdot, mctx->Xdot)); 4319566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 4329566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, mctx->X, INSERT_VALUES, mctx->localX0)); 4339566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, mctx->X, INSERT_VALUES, mctx->localX0)); 434*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 435c4762a1bSJed Brown } 436c4762a1bSJed Brown 437c4762a1bSJed Brown /*TEST 438c4762a1bSJed Brown 439c4762a1bSJed Brown build: 440c4762a1bSJed Brown requires: double !complex adolc 441c4762a1bSJed Brown 442c4762a1bSJed Brown test: 443c4762a1bSJed Brown suffix: 1 444c4762a1bSJed Brown args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian 445c4762a1bSJed Brown output_file: output/adr_ex5adj_mf_1.out 446c4762a1bSJed Brown 447c4762a1bSJed Brown test: 448c4762a1bSJed Brown suffix: 2 449c4762a1bSJed Brown nsize: 4 450c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor 451c4762a1bSJed Brown output_file: output/adr_ex5adj_mf_2.out 452c4762a1bSJed Brown 453c4762a1bSJed Brown TEST*/ 454