xref: /petsc/src/ts/tutorials/ex20opt_ic.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static char help[] = "Solves a ODE-constrained optimization problem -- finding the optimal initial conditions for the van der Pol equation.\n";
2c4762a1bSJed Brown 
30e3d61c9SBarry Smith /*
4c4762a1bSJed Brown   Notes:
5c4762a1bSJed Brown   This code demonstrates how to solve an ODE-constrained optimization problem with TAO, TSAdjoint and TS.
6c4762a1bSJed Brown   The nonlinear problem is written in an ODE equivalent form.
7c4762a1bSJed Brown   The objective is to minimize the difference between observation and model prediction by finding optimal values for initial conditions.
8c4762a1bSJed Brown   The gradient is computed with the discrete adjoint of an implicit method or an explicit method, see ex20adj.c for details.
9c4762a1bSJed Brown */
10c4762a1bSJed Brown 
11c4762a1bSJed Brown #include <petsctao.h>
12c4762a1bSJed Brown #include <petscts.h>
13c4762a1bSJed Brown 
14c4762a1bSJed Brown typedef struct _n_User *User;
15c4762a1bSJed Brown struct _n_User {
16c4762a1bSJed Brown   TS        ts;
17c4762a1bSJed Brown   PetscReal mu;
18c4762a1bSJed Brown   PetscReal next_output;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   /* Sensitivity analysis support */
21c4762a1bSJed Brown   PetscInt  steps;
22c4762a1bSJed Brown   PetscReal ftime;
23c4762a1bSJed Brown   Mat       A;                    /* Jacobian matrix for ODE */
24c4762a1bSJed Brown   Mat       Jacp;                 /* JacobianP matrix for ODE*/
25c4762a1bSJed Brown   Mat       H;                    /* Hessian matrix for optimization */
26c4762a1bSJed Brown   Vec       U, Lambda[1], Mup[1]; /* first-order adjoint variables */
27c4762a1bSJed Brown   Vec       Lambda2[2];           /* second-order adjoint variables */
28c4762a1bSJed Brown   Vec       Ihp1[1];              /* working space for Hessian evaluations */
29c4762a1bSJed Brown   Vec       Dir;                  /* direction vector */
30c4762a1bSJed Brown   PetscReal ob[2];                /* observation used by the cost function */
31c4762a1bSJed Brown   PetscBool implicitform;         /* implicit ODE? */
32c4762a1bSJed Brown };
33c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec, PetscScalar[], User);
34c4762a1bSJed Brown 
35c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE  -------------------- */
36c4762a1bSJed Brown 
379371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx) {
38c4762a1bSJed Brown   User               user = (User)ctx;
39c4762a1bSJed Brown   PetscScalar       *f;
40c4762a1bSJed Brown   const PetscScalar *u;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   PetscFunctionBeginUser;
439566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
449566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
45c4762a1bSJed Brown   f[0] = u[1];
46c4762a1bSJed Brown   f[1] = user->mu * ((1. - u[0] * u[0]) * u[1] - u[0]);
479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
49c4762a1bSJed Brown   PetscFunctionReturn(0);
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
529371c9d4SSatish Balay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx) {
53c4762a1bSJed Brown   User               user     = (User)ctx;
54c4762a1bSJed Brown   PetscReal          mu       = user->mu;
55c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
56c4762a1bSJed Brown   PetscScalar        J[2][2];
57c4762a1bSJed Brown   const PetscScalar *u;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   PetscFunctionBeginUser;
609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
61c4762a1bSJed Brown   J[0][0] = 0;
62c4762a1bSJed Brown   J[1][0] = -mu * (2.0 * u[1] * u[0] + 1.);
63c4762a1bSJed Brown   J[0][1] = 1.0;
64c4762a1bSJed Brown   J[1][1] = mu * (1.0 - u[0] * u[0]);
659566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
68c4762a1bSJed Brown   if (A != B) {
699566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
709566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
71c4762a1bSJed Brown   }
729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
73c4762a1bSJed Brown   PetscFunctionReturn(0);
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
769371c9d4SSatish Balay static PetscErrorCode RHSHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
77c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
78c4762a1bSJed Brown   PetscScalar       *vhv;
79c4762a1bSJed Brown   PetscScalar        dJdU[2][2][2] = {{{0}}};
80c4762a1bSJed Brown   PetscInt           i, j, k;
81c4762a1bSJed Brown   User               user = (User)ctx;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   PetscFunctionBeginUser;
849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
879566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   dJdU[1][0][0] = -2. * user->mu * u[1];
90c4762a1bSJed Brown   dJdU[1][1][0] = -2. * user->mu * u[0];
91c4762a1bSJed Brown   dJdU[1][0][1] = -2. * user->mu * u[0];
92c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
93c4762a1bSJed Brown     vhv[j] = 0;
94c4762a1bSJed Brown     for (k = 0; k < 2; k++)
959371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
96c4762a1bSJed Brown   }
979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE  -------------------- */
105c4762a1bSJed Brown 
1069371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) {
107c4762a1bSJed Brown   User               user = (User)ctx;
108c4762a1bSJed Brown   const PetscScalar *u, *udot;
109c4762a1bSJed Brown   PetscScalar       *f;
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   PetscFunctionBeginUser;
1129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1149566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
115c4762a1bSJed Brown   f[0] = udot[0] - u[1];
116c4762a1bSJed Brown   f[1] = udot[1] - user->mu * ((1.0 - u[0] * u[0]) * u[1] - u[0]);
1179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
120c4762a1bSJed Brown   PetscFunctionReturn(0);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown 
1239371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) {
124c4762a1bSJed Brown   User               user     = (User)ctx;
125c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
126c4762a1bSJed Brown   PetscScalar        J[2][2];
127c4762a1bSJed Brown   const PetscScalar *u;
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   PetscFunctionBeginUser;
1309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1319371c9d4SSatish Balay   J[0][0] = a;
1329371c9d4SSatish Balay   J[0][1] = -1.0;
1339371c9d4SSatish Balay   J[1][0] = user->mu * (1.0 + 2.0 * u[0] * u[1]);
1349371c9d4SSatish Balay   J[1][1] = a - user->mu * (1.0 - u[0] * u[0]);
1359566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
1369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1379566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
139c4762a1bSJed Brown   if (A != B) {
1409566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1419566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
142c4762a1bSJed Brown   }
143c4762a1bSJed Brown   PetscFunctionReturn(0);
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
1479371c9d4SSatish Balay static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) {
148c4762a1bSJed Brown   const PetscScalar *u;
149c4762a1bSJed Brown   PetscReal          tfinal, dt;
150c4762a1bSJed Brown   User               user = (User)ctx;
151c4762a1bSJed Brown   Vec                interpolatedU;
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   PetscFunctionBeginUser;
1549566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
1559566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts, &tfinal));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
1589566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(U, &interpolatedU));
1599566063dSJacob Faibussowitsch     PetscCall(TSInterpolate(ts, user->next_output, interpolatedU));
1609566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(interpolatedU, &u));
1619371c9d4SSatish Balay     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%g] %" PetscInt_FMT " TS %g (dt = %g) X %g %g\n", (double)user->next_output, step, (double)t, (double)dt, (double)PetscRealPart(u[0]), (double)PetscRealPart(u[1])));
1629566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(interpolatedU, &u));
1639566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&interpolatedU));
164c4762a1bSJed Brown     user->next_output += 0.1;
165c4762a1bSJed Brown   }
166c4762a1bSJed Brown   PetscFunctionReturn(0);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
1699371c9d4SSatish Balay static PetscErrorCode IHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
170c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
171c4762a1bSJed Brown   PetscScalar       *vhv;
172c4762a1bSJed Brown   PetscScalar        dJdU[2][2][2] = {{{0}}};
173c4762a1bSJed Brown   PetscInt           i, j, k;
174c4762a1bSJed Brown   User               user = (User)ctx;
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   PetscFunctionBeginUser;
1779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
1799566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
1809566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
181c4762a1bSJed Brown   dJdU[1][0][0] = 2. * user->mu * u[1];
182c4762a1bSJed Brown   dJdU[1][1][0] = 2. * user->mu * u[0];
183c4762a1bSJed Brown   dJdU[1][0][1] = 2. * user->mu * u[0];
184c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
185c4762a1bSJed Brown     vhv[j] = 0;
186c4762a1bSJed Brown     for (k = 0; k < 2; k++)
1879371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
188c4762a1bSJed Brown   }
1899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
1919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
193c4762a1bSJed Brown   PetscFunctionReturn(0);
194c4762a1bSJed Brown }
195c4762a1bSJed Brown 
196c4762a1bSJed Brown /* ------------------ User-defined routines for TAO -------------------------- */
197c4762a1bSJed Brown 
1989371c9d4SSatish Balay static PetscErrorCode FormFunctionGradient(Tao tao, Vec IC, PetscReal *f, Vec G, void *ctx) {
199c4762a1bSJed Brown   User               user_ptr = (User)ctx;
200c4762a1bSJed Brown   TS                 ts       = user_ptr->ts;
201c4762a1bSJed Brown   const PetscScalar *x_ptr;
202c4762a1bSJed Brown   PetscScalar       *y_ptr;
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   PetscFunctionBeginUser;
2059566063dSJacob Faibussowitsch   PetscCall(VecCopy(IC, user_ptr->U)); /* set up the initial condition */
206c4762a1bSJed Brown 
2079566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
2089566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
2099566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.001)); /* can be overwritten by command line options */
2109566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
2119566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, user_ptr->U));
212c4762a1bSJed Brown 
2139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user_ptr->U, &x_ptr));
2149566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Lambda[0], &y_ptr));
215c4762a1bSJed Brown   *f       = (x_ptr[0] - user_ptr->ob[0]) * (x_ptr[0] - user_ptr->ob[0]) + (x_ptr[1] - user_ptr->ob[1]) * (x_ptr[1] - user_ptr->ob[1]);
216c4762a1bSJed Brown   y_ptr[0] = 2. * (x_ptr[0] - user_ptr->ob[0]);
217c4762a1bSJed Brown   y_ptr[1] = 2. * (x_ptr[1] - user_ptr->ob[1]);
2189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Lambda[0], &y_ptr));
2199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user_ptr->U, &x_ptr));
220c4762a1bSJed Brown 
2219566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, user_ptr->Lambda, NULL));
2229566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
2239566063dSJacob Faibussowitsch   PetscCall(VecCopy(user_ptr->Lambda[0], G));
224c4762a1bSJed Brown   PetscFunctionReturn(0);
225c4762a1bSJed Brown }
226c4762a1bSJed Brown 
2279371c9d4SSatish Balay static PetscErrorCode FormHessian(Tao tao, Vec U, Mat H, Mat Hpre, void *ctx) {
228c4762a1bSJed Brown   User           user_ptr = (User)ctx;
229c4762a1bSJed Brown   PetscScalar    harr[2];
230c4762a1bSJed Brown   PetscScalar   *x_ptr;
231c4762a1bSJed Brown   const PetscInt rows[2] = {0, 1};
232c4762a1bSJed Brown   PetscInt       col;
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   PetscFunctionBeginUser;
2359566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, user_ptr->U));
2369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Dir, &x_ptr));
237c4762a1bSJed Brown   x_ptr[0] = 1.;
238c4762a1bSJed Brown   x_ptr[1] = 0.;
2399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Dir, &x_ptr));
2409566063dSJacob Faibussowitsch   PetscCall(Adjoint2(user_ptr->U, harr, user_ptr));
241c4762a1bSJed Brown   col = 0;
2429566063dSJacob Faibussowitsch   PetscCall(MatSetValues(H, 2, rows, 1, &col, harr, INSERT_VALUES));
243c4762a1bSJed Brown 
2449566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, user_ptr->U));
2459566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Dir, &x_ptr));
246c4762a1bSJed Brown   x_ptr[0] = 0.;
247c4762a1bSJed Brown   x_ptr[1] = 1.;
2489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Dir, &x_ptr));
2499566063dSJacob Faibussowitsch   PetscCall(Adjoint2(user_ptr->U, harr, user_ptr));
250c4762a1bSJed Brown   col = 1;
2519566063dSJacob Faibussowitsch   PetscCall(MatSetValues(H, 2, rows, 1, &col, harr, INSERT_VALUES));
252c4762a1bSJed Brown 
2539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
2549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
255c4762a1bSJed Brown   if (H != Hpre) {
2569566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Hpre, MAT_FINAL_ASSEMBLY));
2579566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Hpre, MAT_FINAL_ASSEMBLY));
258c4762a1bSJed Brown   }
259c4762a1bSJed Brown   PetscFunctionReturn(0);
260c4762a1bSJed Brown }
261c4762a1bSJed Brown 
2629371c9d4SSatish Balay static PetscErrorCode MatrixFreeHessian(Tao tao, Vec U, Mat H, Mat Hpre, void *ctx) {
263c4762a1bSJed Brown   User user_ptr = (User)ctx;
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   PetscFunctionBeginUser;
2669566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, user_ptr->U));
267c4762a1bSJed Brown   PetscFunctionReturn(0);
268c4762a1bSJed Brown }
269c4762a1bSJed Brown 
270c4762a1bSJed Brown /* ------------ Routines calculating second-order derivatives -------------- */
271c4762a1bSJed Brown 
2720e3d61c9SBarry Smith /*
273c4762a1bSJed Brown   Compute the Hessian-vector product for the cost function using Second-order adjoint
274c4762a1bSJed Brown */
2759371c9d4SSatish Balay PetscErrorCode Adjoint2(Vec U, PetscScalar arr[], User ctx) {
276c4762a1bSJed Brown   TS           ts = ctx->ts;
277c4762a1bSJed Brown   PetscScalar *x_ptr, *y_ptr;
278c4762a1bSJed Brown   Mat          tlmsen;
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   PetscFunctionBeginUser;
2819566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
282c4762a1bSJed Brown 
2839566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
2849566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
2859566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.001));
2869566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
2879566063dSJacob Faibussowitsch   PetscCall(TSSetCostHessianProducts(ts, 1, ctx->Lambda2, NULL, ctx->Dir));
2889566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetForward(ts, NULL));
2899566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
290c4762a1bSJed Brown 
291c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
2929566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &x_ptr));
2939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda[0], &y_ptr));
294c4762a1bSJed Brown   y_ptr[0] = 2. * (x_ptr[0] - ctx->ob[0]);
295c4762a1bSJed Brown   y_ptr[1] = 2. * (x_ptr[1] - ctx->ob[1]);
2969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda[0], &y_ptr));
2979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &x_ptr));
2989566063dSJacob Faibussowitsch   PetscCall(TSForwardGetSensitivities(ts, NULL, &tlmsen));
2999566063dSJacob Faibussowitsch   PetscCall(MatDenseGetColumn(tlmsen, 0, &x_ptr));
3009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda2[0], &y_ptr));
301c4762a1bSJed Brown   y_ptr[0] = 2. * x_ptr[0];
302c4762a1bSJed Brown   y_ptr[1] = 2. * x_ptr[1];
3039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0], &y_ptr));
3049566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreColumn(tlmsen, &x_ptr));
305c4762a1bSJed Brown 
3069566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, ctx->Lambda, NULL));
307c4762a1bSJed Brown   if (ctx->implicitform) {
3089566063dSJacob Faibussowitsch     PetscCall(TSSetIHessianProduct(ts, ctx->Ihp1, IHessianProductUU, NULL, NULL, NULL, NULL, NULL, NULL, ctx));
309c4762a1bSJed Brown   } else {
3109566063dSJacob Faibussowitsch     PetscCall(TSSetRHSHessianProduct(ts, ctx->Ihp1, RHSHessianProductUU, NULL, NULL, NULL, NULL, NULL, NULL, ctx));
311c4762a1bSJed Brown   }
3129566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
313c4762a1bSJed Brown 
3149566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda2[0], &x_ptr));
315c4762a1bSJed Brown   arr[0] = x_ptr[0];
316c4762a1bSJed Brown   arr[1] = x_ptr[1];
3179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0], &x_ptr));
318c4762a1bSJed Brown 
3199566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
3209566063dSJacob Faibussowitsch   PetscCall(TSAdjointResetForward(ts));
321c4762a1bSJed Brown   PetscFunctionReturn(0);
322c4762a1bSJed Brown }
323c4762a1bSJed Brown 
3249371c9d4SSatish Balay PetscErrorCode FiniteDiff(Vec U, PetscScalar arr[], User ctx) {
325c4762a1bSJed Brown   Vec               Up, G, Gp;
326c4762a1bSJed Brown   const PetscScalar eps = PetscRealConstant(1e-7);
327c4762a1bSJed Brown   PetscScalar      *u;
328c4762a1bSJed Brown   Tao               tao = NULL;
329c4762a1bSJed Brown   PetscReal         f;
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   PetscFunctionBeginUser;
3329566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &Up));
3339566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &G));
3349566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &Gp));
335c4762a1bSJed Brown 
3369566063dSJacob Faibussowitsch   PetscCall(FormFunctionGradient(tao, U, &f, G, ctx));
337c4762a1bSJed Brown 
3389566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, Up));
3399566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Up, &u));
340c4762a1bSJed Brown   u[0] += eps;
3419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Up, &u));
3429566063dSJacob Faibussowitsch   PetscCall(FormFunctionGradient(tao, Up, &f, Gp, ctx));
3439566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Gp, -1, G));
3449566063dSJacob Faibussowitsch   PetscCall(VecScale(Gp, 1. / eps));
3459566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Gp, &u));
346c4762a1bSJed Brown   arr[0] = u[0];
347c4762a1bSJed Brown   arr[1] = u[1];
3489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Gp, &u));
349c4762a1bSJed Brown 
3509566063dSJacob Faibussowitsch   PetscCall(VecCopy(U, Up));
3519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Up, &u));
352c4762a1bSJed Brown   u[1] += eps;
3539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Up, &u));
3549566063dSJacob Faibussowitsch   PetscCall(FormFunctionGradient(tao, Up, &f, Gp, ctx));
3559566063dSJacob Faibussowitsch   PetscCall(VecAXPY(Gp, -1, G));
3569566063dSJacob Faibussowitsch   PetscCall(VecScale(Gp, 1. / eps));
3579566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Gp, &u));
358c4762a1bSJed Brown   arr[2] = u[0];
359c4762a1bSJed Brown   arr[3] = u[1];
3609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Gp, &u));
361c4762a1bSJed Brown 
3629566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&G));
3639566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Gp));
3649566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Up));
365c4762a1bSJed Brown   PetscFunctionReturn(0);
366c4762a1bSJed Brown }
367c4762a1bSJed Brown 
3689371c9d4SSatish Balay static PetscErrorCode HessianProductMat(Mat mat, Vec svec, Vec y) {
369c4762a1bSJed Brown   User         user_ptr;
370c4762a1bSJed Brown   PetscScalar *y_ptr;
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   PetscFunctionBeginUser;
3739566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(mat, &user_ptr));
3749566063dSJacob Faibussowitsch   PetscCall(VecCopy(svec, user_ptr->Dir));
3759566063dSJacob Faibussowitsch   PetscCall(VecGetArray(y, &y_ptr));
3769566063dSJacob Faibussowitsch   PetscCall(Adjoint2(user_ptr->U, y_ptr, user_ptr));
3779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(y, &y_ptr));
378c4762a1bSJed Brown   PetscFunctionReturn(0);
379c4762a1bSJed Brown }
380c4762a1bSJed Brown 
3819371c9d4SSatish Balay int main(int argc, char **argv) {
382c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE, mf = PETSC_TRUE;
383c4762a1bSJed Brown   PetscInt       mode = 0;
384c4762a1bSJed Brown   PetscMPIInt    size;
385c4762a1bSJed Brown   struct _n_User user;
386c4762a1bSJed Brown   Vec            x; /* working vector for TAO */
387c4762a1bSJed Brown   PetscScalar   *x_ptr, arr[4];
388c4762a1bSJed Brown   PetscScalar    ic1 = 2.2, ic2 = -0.7; /* initial guess for TAO */
389c4762a1bSJed Brown   Tao            tao;
390c4762a1bSJed Brown   KSP            ksp;
391c4762a1bSJed Brown   PC             pc;
392c4762a1bSJed Brown 
393c4762a1bSJed Brown   /* Initialize program */
394327415f7SBarry Smith   PetscFunctionBeginUser;
3959566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
3973c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
398c4762a1bSJed Brown 
399c4762a1bSJed Brown   /* Set runtime options */
400c4762a1bSJed Brown   user.next_output  = 0.0;
401c4762a1bSJed Brown   user.mu           = 1.0e3;
402c4762a1bSJed Brown   user.steps        = 0;
403c4762a1bSJed Brown   user.ftime        = 0.5;
404c4762a1bSJed Brown   user.implicitform = PETSC_TRUE;
405c4762a1bSJed Brown 
4069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
4079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
4089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mode", &mode, NULL));
4099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-ic1", &ic1, NULL));
4109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-ic2", &ic2, NULL));
4119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-my_tao_mf", &mf, NULL)); /* matrix-free hessian for optimization */
4129566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &user.implicitform, NULL));
413c4762a1bSJed Brown 
414c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
4159566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.A));
4169566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
4179566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.A));
4189566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.A));
4199566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.U, NULL));
4209566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Dir, NULL));
4219566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Lambda[0], NULL));
4229566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Lambda2[0], NULL));
4239566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Ihp1[0], NULL));
424c4762a1bSJed Brown 
425c4762a1bSJed Brown   /* Create timestepping solver context */
4269566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &user.ts));
4279566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(user.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
428c4762a1bSJed Brown   if (user.implicitform) {
4299566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(user.ts, NULL, IFunction, &user));
4309566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(user.ts, user.A, user.A, IJacobian, &user));
4319566063dSJacob Faibussowitsch     PetscCall(TSSetType(user.ts, TSCN));
432c4762a1bSJed Brown   } else {
4339566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(user.ts, NULL, RHSFunction, &user));
4349566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(user.ts, user.A, user.A, RHSJacobian, &user));
4359566063dSJacob Faibussowitsch     PetscCall(TSSetType(user.ts, TSRK));
436c4762a1bSJed Brown   }
4379566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(user.ts, user.ftime));
4389566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(user.ts, TS_EXACTFINALTIME_MATCHSTEP));
439c4762a1bSJed Brown 
440*48a46eb9SPierre Jolivet   if (monitor) PetscCall(TSMonitorSet(user.ts, Monitor, &user, NULL));
441c4762a1bSJed Brown 
442c4762a1bSJed Brown   /* Set ODE initial conditions */
4439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.U, &x_ptr));
444c4762a1bSJed Brown   x_ptr[0] = 2.0;
445c4762a1bSJed Brown   x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu);
4469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.U, &x_ptr));
447c4762a1bSJed Brown 
448c4762a1bSJed Brown   /* Set runtime options */
4499566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(user.ts));
450c4762a1bSJed Brown 
451c4762a1bSJed Brown   /* Obtain the observation */
4529566063dSJacob Faibussowitsch   PetscCall(TSSolve(user.ts, user.U));
4539566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.U, &x_ptr));
454c4762a1bSJed Brown   user.ob[0] = x_ptr[0];
455c4762a1bSJed Brown   user.ob[1] = x_ptr[1];
4569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.U, &x_ptr));
457c4762a1bSJed Brown 
4589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.U, &x));
4599566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &x_ptr));
460c4762a1bSJed Brown   x_ptr[0] = ic1;
461c4762a1bSJed Brown   x_ptr[1] = ic2;
4629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &x_ptr));
463c4762a1bSJed Brown 
464c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used */
4659566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(user.ts));
466c4762a1bSJed Brown 
467c4762a1bSJed Brown   /* Compare finite difference and second-order adjoint. */
468c4762a1bSJed Brown   switch (mode) {
469c4762a1bSJed Brown   case 2:
4709566063dSJacob Faibussowitsch     PetscCall(FiniteDiff(x, arr, &user));
4719566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n Finite difference approximation of the Hessian\n"));
4729566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%g %g\n%g %g\n", (double)arr[0], (double)arr[1], (double)arr[2], (double)arr[3]));
473c4762a1bSJed Brown     break;
474c4762a1bSJed Brown   case 1: /* Compute the Hessian column by column */
4759566063dSJacob Faibussowitsch     PetscCall(VecCopy(x, user.U));
4769566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user.Dir, &x_ptr));
477c4762a1bSJed Brown     x_ptr[0] = 1.;
478c4762a1bSJed Brown     x_ptr[1] = 0.;
4799566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user.Dir, &x_ptr));
4809566063dSJacob Faibussowitsch     PetscCall(Adjoint2(user.U, arr, &user));
4819566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nFirst column of the Hessian\n"));
4829566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%g\n%g\n", (double)arr[0], (double)arr[1]));
4839566063dSJacob Faibussowitsch     PetscCall(VecCopy(x, user.U));
4849566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user.Dir, &x_ptr));
485c4762a1bSJed Brown     x_ptr[0] = 0.;
486c4762a1bSJed Brown     x_ptr[1] = 1.;
4879566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user.Dir, &x_ptr));
4889566063dSJacob Faibussowitsch     PetscCall(Adjoint2(user.U, arr, &user));
4899566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nSecond column of the Hessian\n"));
4909566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%g\n%g\n", (double)arr[0], (double)arr[1]));
491c4762a1bSJed Brown     break;
492c4762a1bSJed Brown   case 0:
493c4762a1bSJed Brown     /* Create optimization context and set up */
4949566063dSJacob Faibussowitsch     PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
4959566063dSJacob Faibussowitsch     PetscCall(TaoSetType(tao, TAOBLMVM));
4969566063dSJacob Faibussowitsch     PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
497c4762a1bSJed Brown 
498c4762a1bSJed Brown     if (mf) {
4999566063dSJacob Faibussowitsch       PetscCall(MatCreateShell(PETSC_COMM_SELF, 2, 2, 2, 2, (void *)&user, &user.H));
5009566063dSJacob Faibussowitsch       PetscCall(MatShellSetOperation(user.H, MATOP_MULT, (void (*)(void))HessianProductMat));
5019566063dSJacob Faibussowitsch       PetscCall(MatSetOption(user.H, MAT_SYMMETRIC, PETSC_TRUE));
5029566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao, user.H, user.H, MatrixFreeHessian, (void *)&user));
503c4762a1bSJed Brown     } else { /* Create Hessian matrix */
5049566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_WORLD, &user.H));
5059566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(user.H, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
5069566063dSJacob Faibussowitsch       PetscCall(MatSetUp(user.H));
5079566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user));
508c4762a1bSJed Brown     }
509c4762a1bSJed Brown 
510c4762a1bSJed Brown     /* Not use any preconditioner */
5119566063dSJacob Faibussowitsch     PetscCall(TaoGetKSP(tao, &ksp));
512c4762a1bSJed Brown     if (ksp) {
5139566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp, &pc));
5149566063dSJacob Faibussowitsch       PetscCall(PCSetType(pc, PCNONE));
515c4762a1bSJed Brown     }
516c4762a1bSJed Brown 
517c4762a1bSJed Brown     /* Set initial solution guess */
5189566063dSJacob Faibussowitsch     PetscCall(TaoSetSolution(tao, x));
5199566063dSJacob Faibussowitsch     PetscCall(TaoSetFromOptions(tao));
5209566063dSJacob Faibussowitsch     PetscCall(TaoSolve(tao));
5219566063dSJacob Faibussowitsch     PetscCall(TaoDestroy(&tao));
5229566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&user.H));
523c4762a1bSJed Brown     break;
5249371c9d4SSatish Balay   default: break;
525c4762a1bSJed Brown   }
526c4762a1bSJed Brown 
527c4762a1bSJed Brown   /* Free work space.  All PETSc objects should be destroyed when they are no longer needed. */
5289566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.A));
5299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.U));
5309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Lambda[0]));
5319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Lambda2[0]));
5329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp1[0]));
5339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Dir));
5349566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&user.ts));
5359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
5369566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
537b122ec5aSJacob Faibussowitsch   return 0;
538c4762a1bSJed Brown }
539c4762a1bSJed Brown 
540c4762a1bSJed Brown /*TEST
541c4762a1bSJed Brown     build:
542c4762a1bSJed Brown       requires: !complex !single
543c4762a1bSJed Brown 
544c4762a1bSJed Brown     test:
545c4762a1bSJed Brown       args:  -ts_type cn -viewer_binary_skip_info -tao_monitor -tao_view -mu 1000 -ts_dt 0.03125
546c4762a1bSJed Brown       output_file: output/ex20opt_ic_1.out
547c4762a1bSJed Brown 
548c4762a1bSJed Brown     test:
549c4762a1bSJed Brown       suffix: 2
550c4762a1bSJed Brown       args:  -ts_type beuler -viewer_binary_skip_info -tao_monitor -tao_view -mu 100 -ts_dt 0.01 -tao_type bntr -tao_bnk_pc_type none
551c4762a1bSJed Brown       output_file: output/ex20opt_ic_2.out
552c4762a1bSJed Brown 
553c4762a1bSJed Brown     test:
554c4762a1bSJed Brown       suffix: 3
555c4762a1bSJed Brown       args:  -ts_type cn -viewer_binary_skip_info -tao_monitor -tao_view -mu 100 -ts_dt 0.01 -tao_type bntr -tao_bnk_pc_type none
556c4762a1bSJed Brown       output_file: output/ex20opt_ic_3.out
557c4762a1bSJed Brown 
558c4762a1bSJed Brown     test:
559c4762a1bSJed Brown       suffix: 4
560c4762a1bSJed Brown       args: -implicitform 0 -ts_dt 0.01 -ts_max_steps 2 -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view -mode 1 -my_tao_mf
561c4762a1bSJed Brown TEST*/
562