1c4762a1bSJed Brown static char help[] = "Demonstrates adjoint sensitivity analysis for Reaction-Diffusion Equations.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown See ex5.c for details on the equation. 5c4762a1bSJed Brown This code demonestrates the TSAdjoint/TAO interface to solve an inverse initial value problem built on a system of 6c4762a1bSJed Brown time-dependent partial differential equations. 7c4762a1bSJed Brown In this problem, the initial value for the PDE is unknown, but the output (the final solution of the PDE) is known. 8c4762a1bSJed Brown We want to determine the initial value that can produce the given output. 9c4762a1bSJed Brown We formulate the problem as a nonlinear optimization problem that minimizes the discrepency between the simulated 10c4762a1bSJed Brown result and given reference solution, calculate the gradient of the objective function with the discrete adjoint 11c4762a1bSJed Brown solver, and solve the optimization problem with TAO. 12c4762a1bSJed Brown 13c4762a1bSJed Brown Runtime options: 14c4762a1bSJed Brown -forwardonly - run only the forward simulation 15c4762a1bSJed Brown -implicitform - provide IFunction and IJacobian to TS, if not set, RHSFunction and RHSJacobian will be used 16c4762a1bSJed Brown */ 17c4762a1bSJed Brown 1860f0b76eSHong Zhang #include "reaction_diffusion.h" 19c4762a1bSJed Brown #include <petscdm.h> 20c4762a1bSJed Brown #include <petscdmda.h> 21c4762a1bSJed Brown #include <petsctao.h> 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* User-defined routines */ 24c4762a1bSJed Brown extern PetscErrorCode FormFunctionAndGradient(Tao, Vec, PetscReal *, Vec, void *); 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* 27c4762a1bSJed Brown Set terminal condition for the adjoint variable 28c4762a1bSJed Brown */ 299371c9d4SSatish Balay PetscErrorCode InitializeLambda(DM da, Vec lambda, Vec U, AppCtx *appctx) { 30c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = ""; 31c4762a1bSJed Brown PetscViewer viewer; 32c4762a1bSJed Brown Vec Uob; 33c4762a1bSJed Brown 34362febeeSStefano Zampini PetscFunctionBegin; 359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &Uob)); 369566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(filename, sizeof filename, "ex5opt.ob")); 379566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &viewer)); 389566063dSJacob Faibussowitsch PetscCall(VecLoad(Uob, viewer)); 399566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 409566063dSJacob Faibussowitsch PetscCall(VecAYPX(Uob, -1., U)); 419566063dSJacob Faibussowitsch PetscCall(VecScale(Uob, 2.0)); 429566063dSJacob Faibussowitsch PetscCall(VecAXPY(lambda, 1., Uob)); 439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Uob)); 44c4762a1bSJed Brown PetscFunctionReturn(0); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* 48c4762a1bSJed Brown Set up a viewer that dumps data into a binary file 49c4762a1bSJed Brown */ 509371c9d4SSatish Balay PetscErrorCode OutputBIN(DM da, const char *filename, PetscViewer *viewer) { 51362febeeSStefano Zampini PetscFunctionBegin; 529566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)da), viewer)); 539566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERBINARY)); 549566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*viewer, FILE_MODE_WRITE)); 559566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*viewer, filename)); 56c4762a1bSJed Brown PetscFunctionReturn(0); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* 60c4762a1bSJed Brown Generate a reference solution and save it to a binary file 61c4762a1bSJed Brown */ 629371c9d4SSatish Balay PetscErrorCode GenerateOBs(TS ts, Vec U, AppCtx *appctx) { 63c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = ""; 64c4762a1bSJed Brown PetscViewer viewer; 65c4762a1bSJed Brown DM da; 66c4762a1bSJed Brown 67362febeeSStefano Zampini PetscFunctionBegin; 689566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 699566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 709566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(filename, sizeof filename, "ex5opt.ob")); 719566063dSJacob Faibussowitsch PetscCall(OutputBIN(da, filename, &viewer)); 729566063dSJacob Faibussowitsch PetscCall(VecView(U, viewer)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 74c4762a1bSJed Brown PetscFunctionReturn(0); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown 779371c9d4SSatish Balay PetscErrorCode InitialConditions(DM da, Vec U) { 78c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My; 79c4762a1bSJed Brown Field **u; 80c4762a1bSJed Brown PetscReal hx, hy, x, y; 81c4762a1bSJed Brown 82c4762a1bSJed Brown PetscFunctionBegin; 839566063dSJacob 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)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx; 86c4762a1bSJed Brown hy = 2.5 / (PetscReal)My; 87c4762a1bSJed Brown 889566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 89c4762a1bSJed Brown /* Get local grid boundaries */ 909566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 93c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 94c4762a1bSJed Brown y = j * hy; 95c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 96c4762a1bSJed Brown x = i * hx; 97c4762a1bSJed Brown if ((1.0 <= x) && (x <= 1.5) && (1.0 <= y) && (y <= 1.5)) u[j][i].v = .25 * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * x), 2.0) * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * y), 2.0); 98c4762a1bSJed Brown else u[j][i].v = 0.0; 99c4762a1bSJed Brown 100c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v; 101c4762a1bSJed Brown } 102c4762a1bSJed Brown } 103c4762a1bSJed Brown 1049566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 105c4762a1bSJed Brown PetscFunctionReturn(0); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown 1089371c9d4SSatish Balay PetscErrorCode PerturbedInitialConditions(DM da, Vec U) { 109c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My; 110c4762a1bSJed Brown Field **u; 111c4762a1bSJed Brown PetscReal hx, hy, x, y; 112c4762a1bSJed Brown 113c4762a1bSJed Brown PetscFunctionBegin; 1149566063dSJacob 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)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx; 117c4762a1bSJed Brown hy = 2.5 / (PetscReal)My; 118c4762a1bSJed Brown 1199566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 120c4762a1bSJed Brown /* Get local grid boundaries */ 1219566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 124c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 125c4762a1bSJed Brown y = j * hy; 126c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 127c4762a1bSJed Brown x = i * hx; 128c4762a1bSJed Brown if ((1.0 <= x) && (x <= 1.5) && (1.0 <= y) && (y <= 1.5)) u[j][i].v = .25 * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * 0.5 * x), 2.0) * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * 0.5 * y), 2.0); 129c4762a1bSJed Brown else u[j][i].v = 0.0; 130c4762a1bSJed Brown 131c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v; 132c4762a1bSJed Brown } 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 1359566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 136c4762a1bSJed Brown PetscFunctionReturn(0); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 1399371c9d4SSatish Balay PetscErrorCode PerturbedInitialConditions2(DM da, Vec U) { 140c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My; 141c4762a1bSJed Brown Field **u; 142c4762a1bSJed Brown PetscReal hx, hy, x, y; 143c4762a1bSJed Brown 144c4762a1bSJed Brown PetscFunctionBegin; 1459566063dSJacob 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)); 146c4762a1bSJed Brown 147c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx; 148c4762a1bSJed Brown hy = 2.5 / (PetscReal)My; 149c4762a1bSJed Brown 1509566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 151c4762a1bSJed Brown /* Get local grid boundaries */ 1529566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 155c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 156c4762a1bSJed Brown y = j * hy; 157c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 158c4762a1bSJed Brown x = i * hx; 1599371c9d4SSatish Balay if (PetscApproximateGTE(x, 1.0) && PetscApproximateLTE(x, 1.5) && PetscApproximateGTE(y, 1.0) && PetscApproximateLTE(y, 1.5)) 1609371c9d4SSatish Balay u[j][i].v = PetscPowReal(PetscSinReal(4.0 * PETSC_PI * (x - 0.5)), 2.0) * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * y), 2.0) / 4.0; 161c4762a1bSJed Brown else u[j][i].v = 0.0; 162c4762a1bSJed Brown 163c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 1679566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 168c4762a1bSJed Brown PetscFunctionReturn(0); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 1719371c9d4SSatish Balay PetscErrorCode PerturbedInitialConditions3(DM da, Vec U) { 172c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My; 173c4762a1bSJed Brown Field **u; 174c4762a1bSJed Brown PetscReal hx, hy, x, y; 175c4762a1bSJed Brown 176c4762a1bSJed Brown PetscFunctionBegin; 1779566063dSJacob 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)); 178c4762a1bSJed Brown 179c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx; 180c4762a1bSJed Brown hy = 2.5 / (PetscReal)My; 181c4762a1bSJed Brown 1829566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 183c4762a1bSJed Brown /* Get local grid boundaries */ 1849566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 187c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 188c4762a1bSJed Brown y = j * hy; 189c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 190c4762a1bSJed Brown x = i * hx; 191c4762a1bSJed Brown if ((0.5 <= x) && (x <= 2.0) && (0.5 <= y) && (y <= 2.0)) u[j][i].v = .25 * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * x), 2.0) * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * y), 2.0); 192c4762a1bSJed Brown else u[j][i].v = 0.0; 193c4762a1bSJed Brown 194c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown } 197c4762a1bSJed Brown 1989566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 199c4762a1bSJed Brown PetscFunctionReturn(0); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown 2029371c9d4SSatish Balay int main(int argc, char **argv) { 203c4762a1bSJed Brown DM da; 204c4762a1bSJed Brown AppCtx appctx; 205c4762a1bSJed Brown PetscBool forwardonly = PETSC_FALSE, implicitform = PETSC_FALSE; 206c4762a1bSJed Brown PetscInt perturbic = 1; 207c4762a1bSJed Brown 208327415f7SBarry Smith PetscFunctionBeginUser; 2099566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 2109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-forwardonly", &forwardonly, NULL)); 2119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &implicitform, NULL)); 2129566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-perturbic", &perturbic, NULL)); 213c4762a1bSJed Brown 214c4762a1bSJed Brown appctx.D1 = 8.0e-5; 215c4762a1bSJed Brown appctx.D2 = 4.0e-5; 216c4762a1bSJed Brown appctx.gamma = .024; 217c4762a1bSJed Brown appctx.kappa = .06; 21860f0b76eSHong Zhang appctx.aijpc = PETSC_FALSE; 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* Create distributed array (DMDA) to manage parallel grid and vectors */ 2219566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_PERIODIC, DMDA_STENCIL_STAR, 64, 64, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, NULL, &da)); 2229566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 2239566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 2249566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "u")); 2259566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "v")); 226c4762a1bSJed Brown 227c4762a1bSJed Brown /* Extract global vectors from DMDA; then duplicate for remaining 228c4762a1bSJed Brown vectors that are the same types */ 2299566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &appctx.U)); 230c4762a1bSJed Brown 231c4762a1bSJed Brown /* Create timestepping solver context */ 2329566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &appctx.ts)); 2339566063dSJacob Faibussowitsch PetscCall(TSSetType(appctx.ts, TSCN)); 2349566063dSJacob Faibussowitsch PetscCall(TSSetDM(appctx.ts, da)); 2359566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(appctx.ts, TS_NONLINEAR)); 2369566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(appctx.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 237c4762a1bSJed Brown if (!implicitform) { 2389566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(appctx.ts, NULL, RHSFunction, &appctx)); 2399566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(appctx.ts, NULL, NULL, RHSJacobian, &appctx)); 240c4762a1bSJed Brown } else { 2419566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(appctx.ts, NULL, IFunction, &appctx)); 2429566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(appctx.ts, NULL, NULL, IJacobian, &appctx)); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown 245c4762a1bSJed Brown /* Set initial conditions */ 2469566063dSJacob Faibussowitsch PetscCall(InitialConditions(da, appctx.U)); 2479566063dSJacob Faibussowitsch PetscCall(TSSetSolution(appctx.ts, appctx.U)); 248c4762a1bSJed Brown 249c4762a1bSJed Brown /* Set solver options */ 2509566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(appctx.ts, 2000.0)); 2519566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(appctx.ts, 0.5)); 2529566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(appctx.ts, TS_EXACTFINALTIME_MATCHSTEP)); 2539566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(appctx.ts)); 254c4762a1bSJed Brown 2559566063dSJacob Faibussowitsch PetscCall(GenerateOBs(appctx.ts, appctx.U, &appctx)); 256c4762a1bSJed Brown 257c4762a1bSJed Brown if (!forwardonly) { 258c4762a1bSJed Brown Tao tao; 259c4762a1bSJed Brown Vec P; 260c4762a1bSJed Brown Vec lambda[1]; 261956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 262c4762a1bSJed Brown PetscLogStage opt_stage; 263956f8c0dSBarry Smith #endif 264c4762a1bSJed Brown 2659566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Optimization", &opt_stage)); 2669566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(opt_stage)); 267c4762a1bSJed Brown if (perturbic == 1) { 2689566063dSJacob Faibussowitsch PetscCall(PerturbedInitialConditions(da, appctx.U)); 269c4762a1bSJed Brown } else if (perturbic == 2) { 2709566063dSJacob Faibussowitsch PetscCall(PerturbedInitialConditions2(da, appctx.U)); 271c4762a1bSJed Brown } else if (perturbic == 3) { 2729566063dSJacob Faibussowitsch PetscCall(PerturbedInitialConditions3(da, appctx.U)); 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 2759566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx.U, &lambda[0])); 2769566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(appctx.ts, 1, lambda, NULL)); 277c4762a1bSJed Brown 278c4762a1bSJed Brown /* Have the TS save its trajectory needed by TSAdjointSolve() */ 2799566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(appctx.ts)); 280c4762a1bSJed Brown 281c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 2829566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 2839566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBLMVM)); 284c4762a1bSJed Brown 285c4762a1bSJed Brown /* Set initial guess for TAO */ 2869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx.U, &P)); 2879566063dSJacob Faibussowitsch PetscCall(VecCopy(appctx.U, P)); 2889566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, P)); 289c4762a1bSJed Brown 290c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 2919566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionAndGradient, &appctx)); 292c4762a1bSJed Brown 293c4762a1bSJed Brown /* Check for any TAO command line options */ 2949566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 295c4762a1bSJed Brown 2969566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 2979566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 2989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0])); 2999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&P)); 3009566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 301c4762a1bSJed Brown } 302c4762a1bSJed Brown 303c4762a1bSJed Brown /* Free work space. All PETSc objects should be destroyed when they 304c4762a1bSJed Brown are no longer needed. */ 3059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.U)); 3069566063dSJacob Faibussowitsch PetscCall(TSDestroy(&appctx.ts)); 3079566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 3089566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 309b122ec5aSJacob Faibussowitsch return 0; 310c4762a1bSJed Brown } 311c4762a1bSJed Brown 312c4762a1bSJed Brown /* ------------------------ TAO callbacks ---------------------------- */ 313c4762a1bSJed Brown 314c4762a1bSJed Brown /* 315c4762a1bSJed Brown FormFunctionAndGradient - Evaluates the function and corresponding gradient. 316c4762a1bSJed Brown 317c4762a1bSJed Brown Input Parameters: 318c4762a1bSJed Brown tao - the Tao context 319c4762a1bSJed Brown P - the input vector 320a82e8c82SStefano Zampini ctx - optional user-defined context, as set by TaoSetObjectiveAndGradient() 321c4762a1bSJed Brown 322c4762a1bSJed Brown Output Parameters: 323c4762a1bSJed Brown f - the newly evaluated function 324c4762a1bSJed Brown G - the newly evaluated gradient 325c4762a1bSJed Brown */ 3269371c9d4SSatish Balay PetscErrorCode FormFunctionAndGradient(Tao tao, Vec P, PetscReal *f, Vec G, void *ctx) { 327c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; 328c4762a1bSJed Brown PetscReal soberr, timestep; 329c4762a1bSJed Brown Vec *lambda; 330c4762a1bSJed Brown Vec SDiff; 331c4762a1bSJed Brown DM da; 332c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = ""; 333c4762a1bSJed Brown PetscViewer viewer; 334c4762a1bSJed Brown 335c4762a1bSJed Brown PetscFunctionBeginUser; 3369566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx->ts, 0.0)); 3379566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(appctx->ts, ×tep)); 338*48a46eb9SPierre Jolivet if (timestep < 0) PetscCall(TSSetTimeStep(appctx->ts, -timestep)); 3399566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(appctx->ts, 0)); 3409566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(appctx->ts)); 341c4762a1bSJed Brown 3429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(P, &SDiff)); 3439566063dSJacob Faibussowitsch PetscCall(VecCopy(P, appctx->U)); 3449566063dSJacob Faibussowitsch PetscCall(TSGetDM(appctx->ts, &da)); 345c4762a1bSJed Brown *f = 0; 346c4762a1bSJed Brown 3479566063dSJacob Faibussowitsch PetscCall(TSSolve(appctx->ts, appctx->U)); 3489566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(filename, sizeof filename, "ex5opt.ob")); 3499566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &viewer)); 3509566063dSJacob Faibussowitsch PetscCall(VecLoad(SDiff, viewer)); 3519566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 3529566063dSJacob Faibussowitsch PetscCall(VecAYPX(SDiff, -1., appctx->U)); 3539566063dSJacob Faibussowitsch PetscCall(VecDot(SDiff, SDiff, &soberr)); 354c4762a1bSJed Brown *f += soberr; 355c4762a1bSJed Brown 3569566063dSJacob Faibussowitsch PetscCall(TSGetCostGradients(appctx->ts, NULL, &lambda, NULL)); 3579566063dSJacob Faibussowitsch PetscCall(VecSet(lambda[0], 0.0)); 3589566063dSJacob Faibussowitsch PetscCall(InitializeLambda(da, lambda[0], appctx->U, appctx)); 3599566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(appctx->ts)); 360c4762a1bSJed Brown 3619566063dSJacob Faibussowitsch PetscCall(VecCopy(lambda[0], G)); 362c4762a1bSJed Brown 3639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&SDiff)); 364c4762a1bSJed Brown PetscFunctionReturn(0); 365c4762a1bSJed Brown } 366c4762a1bSJed Brown 367c4762a1bSJed Brown /*TEST 368c4762a1bSJed Brown 369c4762a1bSJed Brown build: 37060f0b76eSHong Zhang depends: reaction_diffusion.c 371c4762a1bSJed Brown requires: !complex !single 372c4762a1bSJed Brown 373c4762a1bSJed Brown test: 374c4762a1bSJed Brown args: -ts_max_steps 5 -ts_type rk -ts_rk_type 3 -ts_trajectory_type memory -tao_monitor -tao_view -tao_gatol 1e-6 375c4762a1bSJed Brown output_file: output/ex5opt_ic_1.out 376c4762a1bSJed Brown 377c4762a1bSJed Brown TEST*/ 378