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.
9da81f932SPierre Jolivet We formulate the problem as a nonlinear optimization problem that minimizes the discrepancy 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 */
InitializeLambda(DM da,Vec lambda,Vec U,AppCtx * appctx)29d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeLambda(DM da, Vec lambda, Vec U, AppCtx *appctx)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "";
32c4762a1bSJed Brown PetscViewer viewer;
33c4762a1bSJed Brown Vec Uob;
34c4762a1bSJed Brown
35362febeeSStefano Zampini PetscFunctionBegin;
369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U, &Uob));
379566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(filename, sizeof filename, "ex5opt.ob"));
389566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &viewer));
399566063dSJacob Faibussowitsch PetscCall(VecLoad(Uob, viewer));
409566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
419566063dSJacob Faibussowitsch PetscCall(VecAYPX(Uob, -1., U));
429566063dSJacob Faibussowitsch PetscCall(VecScale(Uob, 2.0));
439566063dSJacob Faibussowitsch PetscCall(VecAXPY(lambda, 1., Uob));
449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Uob));
453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
46c4762a1bSJed Brown }
47c4762a1bSJed Brown
48c4762a1bSJed Brown /*
49c4762a1bSJed Brown Set up a viewer that dumps data into a binary file
50c4762a1bSJed Brown */
OutputBIN(DM da,const char * filename,PetscViewer * viewer)51d71ae5a4SJacob Faibussowitsch PetscErrorCode OutputBIN(DM da, const char *filename, PetscViewer *viewer)
52d71ae5a4SJacob Faibussowitsch {
53362febeeSStefano Zampini PetscFunctionBegin;
549566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)da), viewer));
559566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*viewer, PETSCVIEWERBINARY));
569566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*viewer, FILE_MODE_WRITE));
579566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*viewer, filename));
583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
59c4762a1bSJed Brown }
60c4762a1bSJed Brown
61c4762a1bSJed Brown /*
62c4762a1bSJed Brown Generate a reference solution and save it to a binary file
63c4762a1bSJed Brown */
GenerateOBs(TS ts,Vec U,AppCtx * appctx)64d71ae5a4SJacob Faibussowitsch PetscErrorCode GenerateOBs(TS ts, Vec U, AppCtx *appctx)
65d71ae5a4SJacob Faibussowitsch {
66c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "";
67c4762a1bSJed Brown PetscViewer viewer;
68c4762a1bSJed Brown DM da;
69c4762a1bSJed Brown
70362febeeSStefano Zampini PetscFunctionBegin;
719566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da));
729566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U));
739566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(filename, sizeof filename, "ex5opt.ob"));
749566063dSJacob Faibussowitsch PetscCall(OutputBIN(da, filename, &viewer));
759566063dSJacob Faibussowitsch PetscCall(VecView(U, viewer));
769566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
78c4762a1bSJed Brown }
79c4762a1bSJed Brown
InitialConditions(DM da,Vec U)80d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialConditions(DM da, Vec U)
81d71ae5a4SJacob Faibussowitsch {
82c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My;
83c4762a1bSJed Brown Field **u;
84c4762a1bSJed Brown PetscReal hx, hy, x, y;
85c4762a1bSJed Brown
86c4762a1bSJed Brown PetscFunctionBegin;
879566063dSJacob 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));
88c4762a1bSJed Brown
89c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx;
90c4762a1bSJed Brown hy = 2.5 / (PetscReal)My;
91c4762a1bSJed Brown
929566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u));
93c4762a1bSJed Brown /* Get local grid boundaries */
949566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
95c4762a1bSJed Brown
96c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
97c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
98c4762a1bSJed Brown y = j * hy;
99c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
100c4762a1bSJed Brown x = i * hx;
101c4762a1bSJed 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);
102c4762a1bSJed Brown else u[j][i].v = 0.0;
103c4762a1bSJed Brown
104c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v;
105c4762a1bSJed Brown }
106c4762a1bSJed Brown }
107c4762a1bSJed Brown
1089566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u));
1093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
110c4762a1bSJed Brown }
111c4762a1bSJed Brown
PerturbedInitialConditions(DM da,Vec U)112d71ae5a4SJacob Faibussowitsch PetscErrorCode PerturbedInitialConditions(DM da, Vec U)
113d71ae5a4SJacob Faibussowitsch {
114c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My;
115c4762a1bSJed Brown Field **u;
116c4762a1bSJed Brown PetscReal hx, hy, x, y;
117c4762a1bSJed Brown
118c4762a1bSJed Brown PetscFunctionBegin;
1199566063dSJacob 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));
120c4762a1bSJed Brown
121c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx;
122c4762a1bSJed Brown hy = 2.5 / (PetscReal)My;
123c4762a1bSJed Brown
1249566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u));
125c4762a1bSJed Brown /* Get local grid boundaries */
1269566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
127c4762a1bSJed Brown
128c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
129c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
130c4762a1bSJed Brown y = j * hy;
131c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
132c4762a1bSJed Brown x = i * hx;
133c4762a1bSJed 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);
134c4762a1bSJed Brown else u[j][i].v = 0.0;
135c4762a1bSJed Brown
136c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v;
137c4762a1bSJed Brown }
138c4762a1bSJed Brown }
139c4762a1bSJed Brown
1409566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u));
1413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
142c4762a1bSJed Brown }
143c4762a1bSJed Brown
PerturbedInitialConditions2(DM da,Vec U)144d71ae5a4SJacob Faibussowitsch PetscErrorCode PerturbedInitialConditions2(DM da, Vec U)
145d71ae5a4SJacob Faibussowitsch {
146c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My;
147c4762a1bSJed Brown Field **u;
148c4762a1bSJed Brown PetscReal hx, hy, x, y;
149c4762a1bSJed Brown
150c4762a1bSJed Brown PetscFunctionBegin;
1519566063dSJacob 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));
152c4762a1bSJed Brown
153c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx;
154c4762a1bSJed Brown hy = 2.5 / (PetscReal)My;
155c4762a1bSJed Brown
1569566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u));
157c4762a1bSJed Brown /* Get local grid boundaries */
1589566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
159c4762a1bSJed Brown
160c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
161c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
162c4762a1bSJed Brown y = j * hy;
163c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
164c4762a1bSJed Brown x = i * hx;
1659371c9d4SSatish Balay if (PetscApproximateGTE(x, 1.0) && PetscApproximateLTE(x, 1.5) && PetscApproximateGTE(y, 1.0) && PetscApproximateLTE(y, 1.5))
1669371c9d4SSatish 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;
167c4762a1bSJed Brown else u[j][i].v = 0.0;
168c4762a1bSJed Brown
169c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v;
170c4762a1bSJed Brown }
171c4762a1bSJed Brown }
172c4762a1bSJed Brown
1739566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u));
1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
175c4762a1bSJed Brown }
176c4762a1bSJed Brown
PerturbedInitialConditions3(DM da,Vec U)177d71ae5a4SJacob Faibussowitsch PetscErrorCode PerturbedInitialConditions3(DM da, Vec U)
178d71ae5a4SJacob Faibussowitsch {
179c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My;
180c4762a1bSJed Brown Field **u;
181c4762a1bSJed Brown PetscReal hx, hy, x, y;
182c4762a1bSJed Brown
183c4762a1bSJed Brown PetscFunctionBegin;
1849566063dSJacob 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));
185c4762a1bSJed Brown
186c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx;
187c4762a1bSJed Brown hy = 2.5 / (PetscReal)My;
188c4762a1bSJed Brown
1899566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u));
190c4762a1bSJed Brown /* Get local grid boundaries */
1919566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
192c4762a1bSJed Brown
193c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */
194c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) {
195c4762a1bSJed Brown y = j * hy;
196c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) {
197c4762a1bSJed Brown x = i * hx;
198c4762a1bSJed 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);
199c4762a1bSJed Brown else u[j][i].v = 0.0;
200c4762a1bSJed Brown
201c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v;
202c4762a1bSJed Brown }
203c4762a1bSJed Brown }
204c4762a1bSJed Brown
2059566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u));
2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
207c4762a1bSJed Brown }
208c4762a1bSJed Brown
main(int argc,char ** argv)209d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
210d71ae5a4SJacob Faibussowitsch {
211c4762a1bSJed Brown DM da;
212c4762a1bSJed Brown AppCtx appctx;
213c4762a1bSJed Brown PetscBool forwardonly = PETSC_FALSE, implicitform = PETSC_FALSE;
214c4762a1bSJed Brown PetscInt perturbic = 1;
215c4762a1bSJed Brown
216327415f7SBarry Smith PetscFunctionBeginUser;
217c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-forwardonly", &forwardonly, NULL));
2199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &implicitform, NULL));
2209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-perturbic", &perturbic, NULL));
221c4762a1bSJed Brown
222c4762a1bSJed Brown appctx.D1 = 8.0e-5;
223c4762a1bSJed Brown appctx.D2 = 4.0e-5;
224c4762a1bSJed Brown appctx.gamma = .024;
225c4762a1bSJed Brown appctx.kappa = .06;
22660f0b76eSHong Zhang appctx.aijpc = PETSC_FALSE;
227c4762a1bSJed Brown
228c4762a1bSJed Brown /* Create distributed array (DMDA) to manage parallel grid and vectors */
2299566063dSJacob 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));
2309566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da));
2319566063dSJacob Faibussowitsch PetscCall(DMSetUp(da));
2329566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "u"));
2339566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "v"));
234c4762a1bSJed Brown
235c4762a1bSJed Brown /* Extract global vectors from DMDA; then duplicate for remaining
236c4762a1bSJed Brown vectors that are the same types */
2379566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &appctx.U));
238c4762a1bSJed Brown
239c4762a1bSJed Brown /* Create timestepping solver context */
2409566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &appctx.ts));
2419566063dSJacob Faibussowitsch PetscCall(TSSetType(appctx.ts, TSCN));
2429566063dSJacob Faibussowitsch PetscCall(TSSetDM(appctx.ts, da));
2439566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(appctx.ts, TS_NONLINEAR));
2449566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(appctx.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
245c4762a1bSJed Brown if (!implicitform) {
2469566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(appctx.ts, NULL, RHSFunction, &appctx));
2479566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(appctx.ts, NULL, NULL, RHSJacobian, &appctx));
248c4762a1bSJed Brown } else {
2499566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(appctx.ts, NULL, IFunction, &appctx));
2509566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(appctx.ts, NULL, NULL, IJacobian, &appctx));
251c4762a1bSJed Brown }
252c4762a1bSJed Brown
253c4762a1bSJed Brown /* Set initial conditions */
2549566063dSJacob Faibussowitsch PetscCall(InitialConditions(da, appctx.U));
2559566063dSJacob Faibussowitsch PetscCall(TSSetSolution(appctx.ts, appctx.U));
256c4762a1bSJed Brown
257c4762a1bSJed Brown /* Set solver options */
2589566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(appctx.ts, 2000.0));
2599566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(appctx.ts, 0.5));
2609566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(appctx.ts, TS_EXACTFINALTIME_MATCHSTEP));
2619566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(appctx.ts));
262c4762a1bSJed Brown
2639566063dSJacob Faibussowitsch PetscCall(GenerateOBs(appctx.ts, appctx.U, &appctx));
264c4762a1bSJed Brown
265c4762a1bSJed Brown if (!forwardonly) {
266c4762a1bSJed Brown Tao tao;
267c4762a1bSJed Brown Vec P;
268c4762a1bSJed Brown Vec lambda[1];
269c4762a1bSJed Brown PetscLogStage opt_stage;
270c4762a1bSJed Brown
2719566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Optimization", &opt_stage));
2729566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(opt_stage));
273c4762a1bSJed Brown if (perturbic == 1) {
2749566063dSJacob Faibussowitsch PetscCall(PerturbedInitialConditions(da, appctx.U));
275c4762a1bSJed Brown } else if (perturbic == 2) {
2769566063dSJacob Faibussowitsch PetscCall(PerturbedInitialConditions2(da, appctx.U));
277c4762a1bSJed Brown } else if (perturbic == 3) {
2789566063dSJacob Faibussowitsch PetscCall(PerturbedInitialConditions3(da, appctx.U));
279c4762a1bSJed Brown }
280c4762a1bSJed Brown
2819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx.U, &lambda[0]));
2829566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(appctx.ts, 1, lambda, NULL));
283c4762a1bSJed Brown
284c4762a1bSJed Brown /* Have the TS save its trajectory needed by TSAdjointSolve() */
2859566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(appctx.ts));
286c4762a1bSJed Brown
287c4762a1bSJed Brown /* Create TAO solver and set desired solution method */
2889566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
2899566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBLMVM));
290c4762a1bSJed Brown
291c4762a1bSJed Brown /* Set initial guess for TAO */
2929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(appctx.U, &P));
2939566063dSJacob Faibussowitsch PetscCall(VecCopy(appctx.U, P));
2949566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, P));
295c4762a1bSJed Brown
296c4762a1bSJed Brown /* Set routine for function and gradient evaluation */
2979566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionAndGradient, &appctx));
298c4762a1bSJed Brown
299c4762a1bSJed Brown /* Check for any TAO command line options */
3009566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
301c4762a1bSJed Brown
3029566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
3039566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
3049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0]));
3059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&P));
3069566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop());
307c4762a1bSJed Brown }
308c4762a1bSJed Brown
309c4762a1bSJed Brown /* Free work space. All PETSc objects should be destroyed when they
310c4762a1bSJed Brown are no longer needed. */
3119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.U));
3129566063dSJacob Faibussowitsch PetscCall(TSDestroy(&appctx.ts));
3139566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da));
3149566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
315b122ec5aSJacob Faibussowitsch return 0;
316c4762a1bSJed Brown }
317c4762a1bSJed Brown
318c4762a1bSJed Brown /* ------------------------ TAO callbacks ---------------------------- */
319c4762a1bSJed Brown
320c4762a1bSJed Brown /*
321c4762a1bSJed Brown FormFunctionAndGradient - Evaluates the function and corresponding gradient.
322c4762a1bSJed Brown
323c4762a1bSJed Brown Input Parameters:
324c4762a1bSJed Brown tao - the Tao context
325c4762a1bSJed Brown P - the input vector
326a82e8c82SStefano Zampini ctx - optional user-defined context, as set by TaoSetObjectiveAndGradient()
327c4762a1bSJed Brown
328c4762a1bSJed Brown Output Parameters:
329c4762a1bSJed Brown f - the newly evaluated function
330c4762a1bSJed Brown G - the newly evaluated gradient
331c4762a1bSJed Brown */
FormFunctionAndGradient(Tao tao,Vec P,PetscReal * f,Vec G,PetscCtx ctx)332*2a8381b2SBarry Smith PetscErrorCode FormFunctionAndGradient(Tao tao, Vec P, PetscReal *f, Vec G, PetscCtx ctx)
333d71ae5a4SJacob Faibussowitsch {
334c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx;
335c4762a1bSJed Brown PetscReal soberr, timestep;
336c4762a1bSJed Brown Vec *lambda;
337c4762a1bSJed Brown Vec SDiff;
338c4762a1bSJed Brown DM da;
339c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN] = "";
340c4762a1bSJed Brown PetscViewer viewer;
341c4762a1bSJed Brown
342c4762a1bSJed Brown PetscFunctionBeginUser;
3439566063dSJacob Faibussowitsch PetscCall(TSSetTime(appctx->ts, 0.0));
3449566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(appctx->ts, ×tep));
34548a46eb9SPierre Jolivet if (timestep < 0) PetscCall(TSSetTimeStep(appctx->ts, -timestep));
3469566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(appctx->ts, 0));
3479566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(appctx->ts));
348c4762a1bSJed Brown
3499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(P, &SDiff));
3509566063dSJacob Faibussowitsch PetscCall(VecCopy(P, appctx->U));
3519566063dSJacob Faibussowitsch PetscCall(TSGetDM(appctx->ts, &da));
352c4762a1bSJed Brown *f = 0;
353c4762a1bSJed Brown
3549566063dSJacob Faibussowitsch PetscCall(TSSolve(appctx->ts, appctx->U));
3559566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(filename, sizeof filename, "ex5opt.ob"));
3569566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, filename, FILE_MODE_READ, &viewer));
3579566063dSJacob Faibussowitsch PetscCall(VecLoad(SDiff, viewer));
3589566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer));
3599566063dSJacob Faibussowitsch PetscCall(VecAYPX(SDiff, -1., appctx->U));
3609566063dSJacob Faibussowitsch PetscCall(VecDot(SDiff, SDiff, &soberr));
361c4762a1bSJed Brown *f += soberr;
362c4762a1bSJed Brown
3639566063dSJacob Faibussowitsch PetscCall(TSGetCostGradients(appctx->ts, NULL, &lambda, NULL));
3649566063dSJacob Faibussowitsch PetscCall(VecSet(lambda[0], 0.0));
3659566063dSJacob Faibussowitsch PetscCall(InitializeLambda(da, lambda[0], appctx->U, appctx));
3669566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(appctx->ts));
367c4762a1bSJed Brown
3689566063dSJacob Faibussowitsch PetscCall(VecCopy(lambda[0], G));
369c4762a1bSJed Brown
3709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&SDiff));
3713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
372c4762a1bSJed Brown }
373c4762a1bSJed Brown
374c4762a1bSJed Brown /*TEST
375c4762a1bSJed Brown
376c4762a1bSJed Brown build:
37760f0b76eSHong Zhang depends: reaction_diffusion.c
378c4762a1bSJed Brown requires: !complex !single
379c4762a1bSJed Brown
380c4762a1bSJed Brown test:
381c4762a1bSJed 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
382c4762a1bSJed Brown output_file: output/ex5opt_ic_1.out
383c4762a1bSJed Brown
384c4762a1bSJed Brown TEST*/
385