xref: /petsc/src/ts/tutorials/ex20opt_p.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation.\n\
3c4762a1bSJed Brown Input parameters include:\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /* ------------------------------------------------------------------------
6c4762a1bSJed Brown 
7c4762a1bSJed Brown   Notes:
8c4762a1bSJed Brown   This code demonstrates how to solve a DAE-constrained optimization problem with TAO, TSAdjoint and TS.
9c4762a1bSJed Brown   The nonlinear problem is written in a DAE equivalent form.
10c4762a1bSJed Brown   The objective is to minimize the difference between observation and model prediction by finding an optimal value for parameter \mu.
11c4762a1bSJed Brown   The gradient is computed with the discrete adjoint of an implicit theta method, see ex20adj.c for details.
12c4762a1bSJed Brown   ------------------------------------------------------------------------- */
13c4762a1bSJed Brown #include <petsctao.h>
14c4762a1bSJed Brown #include <petscts.h>
15c4762a1bSJed Brown 
16c4762a1bSJed Brown typedef struct _n_User *User;
17c4762a1bSJed Brown struct _n_User {
18c4762a1bSJed Brown   TS        ts;
19c4762a1bSJed Brown   PetscReal mu;
20c4762a1bSJed Brown   PetscReal next_output;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   /* Sensitivity analysis support */
23c4762a1bSJed Brown   PetscReal ftime;
24c4762a1bSJed Brown   Mat       A;                    /* Jacobian matrix */
25c4762a1bSJed Brown   Mat       Jacp;                 /* JacobianP matrix */
26c4762a1bSJed Brown   Mat       H;                    /* Hessian matrix for optimization */
27c4762a1bSJed Brown   Vec       U, Lambda[1], Mup[1]; /* adjoint variables */
28c4762a1bSJed Brown   Vec       Lambda2[1], Mup2[1];  /* second-order adjoint variables */
29c4762a1bSJed Brown   Vec       Ihp1[1];              /* working space for Hessian evaluations */
30c4762a1bSJed Brown   Vec       Ihp2[1];              /* working space for Hessian evaluations */
31c4762a1bSJed Brown   Vec       Ihp3[1];              /* working space for Hessian evaluations */
32c4762a1bSJed Brown   Vec       Ihp4[1];              /* working space for Hessian evaluations */
33c4762a1bSJed Brown   Vec       Dir;                  /* direction vector */
34c4762a1bSJed Brown   PetscReal ob[2];                /* observation used by the cost function */
35c4762a1bSJed Brown   PetscBool implicitform;         /* implicit ODE? */
36c4762a1bSJed Brown };
37c4762a1bSJed Brown 
38c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
39c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
40c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec, PetscScalar[], User);
41c4762a1bSJed Brown 
42c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE  -------------------- */
43c4762a1bSJed Brown 
449371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx) {
45c4762a1bSJed Brown   User               user = (User)ctx;
46c4762a1bSJed Brown   PetscScalar       *f;
47c4762a1bSJed Brown   const PetscScalar *u;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   PetscFunctionBeginUser;
509566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
52c4762a1bSJed Brown   f[0] = u[1];
53c4762a1bSJed Brown   f[1] = user->mu * ((1. - u[0] * u[0]) * u[1] - u[0]);
549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
56c4762a1bSJed Brown   PetscFunctionReturn(0);
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
599371c9d4SSatish Balay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx) {
60c4762a1bSJed Brown   User               user     = (User)ctx;
61c4762a1bSJed Brown   PetscReal          mu       = user->mu;
62c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
63c4762a1bSJed Brown   PetscScalar        J[2][2];
64c4762a1bSJed Brown   const PetscScalar *u;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   PetscFunctionBeginUser;
679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   J[0][0] = 0;
70c4762a1bSJed Brown   J[1][0] = -mu * (2.0 * u[1] * u[0] + 1.);
71c4762a1bSJed Brown   J[0][1] = 1.0;
72c4762a1bSJed Brown   J[1][1] = mu * (1.0 - u[0] * u[0]);
739566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
76c4762a1bSJed Brown   if (B && A != B) {
779566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
789566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
79c4762a1bSJed Brown   }
80c4762a1bSJed Brown 
819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
82c4762a1bSJed Brown   PetscFunctionReturn(0);
83c4762a1bSJed Brown }
84c4762a1bSJed Brown 
859371c9d4SSatish Balay static PetscErrorCode RHSHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
86c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
87c4762a1bSJed Brown   PetscScalar       *vhv;
88c4762a1bSJed Brown   PetscScalar        dJdU[2][2][2] = {{{0}}};
89c4762a1bSJed Brown   PetscInt           i, j, k;
90c4762a1bSJed Brown   User               user = (User)ctx;
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   PetscFunctionBeginUser;
939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
959566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
969566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   dJdU[1][0][0] = -2. * user->mu * u[1];
99c4762a1bSJed Brown   dJdU[1][1][0] = -2. * user->mu * u[0];
100c4762a1bSJed Brown   dJdU[1][0][1] = -2. * user->mu * u[0];
101c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
102c4762a1bSJed Brown     vhv[j] = 0;
103c4762a1bSJed Brown     for (k = 0; k < 2; k++)
1049371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
105c4762a1bSJed Brown   }
106c4762a1bSJed Brown 
1079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
1099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
111c4762a1bSJed Brown   PetscFunctionReturn(0);
112c4762a1bSJed Brown }
113c4762a1bSJed Brown 
1149371c9d4SSatish Balay static PetscErrorCode RHSHessianProductUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
115c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
116c4762a1bSJed Brown   PetscScalar       *vhv;
117c4762a1bSJed Brown   PetscScalar        dJdP[2][2][1] = {{{0}}};
118c4762a1bSJed Brown   PetscInt           i, j, k;
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   PetscFunctionBeginUser;
1219566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1229566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
1239566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
1249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   dJdP[1][0][0] = -(1. + 2. * u[0] * u[1]);
127c4762a1bSJed Brown   dJdP[1][1][0] = 1. - u[0] * u[0];
128c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
129c4762a1bSJed Brown     vhv[j] = 0;
130c4762a1bSJed Brown     for (k = 0; k < 1; k++)
1319371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdP[i][j][k] * vr[k];
132c4762a1bSJed Brown   }
133c4762a1bSJed Brown 
1349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
1369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
138c4762a1bSJed Brown   PetscFunctionReturn(0);
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
1419371c9d4SSatish Balay static PetscErrorCode RHSHessianProductPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
142c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
143c4762a1bSJed Brown   PetscScalar       *vhv;
144c4762a1bSJed Brown   PetscScalar        dJdU[2][1][2] = {{{0}}};
145c4762a1bSJed Brown   PetscInt           i, j, k;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   PetscFunctionBeginUser;
1489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
1509566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
1519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   dJdU[1][0][0] = -1. - 2. * u[1] * u[0];
154c4762a1bSJed Brown   dJdU[1][0][1] = 1. - u[0] * u[0];
155c4762a1bSJed Brown   for (j = 0; j < 1; j++) {
156c4762a1bSJed Brown     vhv[j] = 0;
157c4762a1bSJed Brown     for (k = 0; k < 2; k++)
1589371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
159c4762a1bSJed Brown   }
160c4762a1bSJed Brown 
1619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
1639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
165c4762a1bSJed Brown   PetscFunctionReturn(0);
166c4762a1bSJed Brown }
167c4762a1bSJed Brown 
1689371c9d4SSatish Balay static PetscErrorCode RHSHessianProductPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
169c4762a1bSJed Brown   PetscFunctionBeginUser;
170c4762a1bSJed Brown   PetscFunctionReturn(0);
171c4762a1bSJed Brown }
172c4762a1bSJed Brown 
173c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE  -------------------- */
174c4762a1bSJed Brown 
1759371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) {
176c4762a1bSJed Brown   User               user = (User)ctx;
177c4762a1bSJed Brown   PetscScalar       *f;
178c4762a1bSJed Brown   const PetscScalar *u, *udot;
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   PetscFunctionBeginUser;
1819566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1829566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1839566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   f[0] = udot[0] - u[1];
186c4762a1bSJed Brown   f[1] = udot[1] - user->mu * ((1.0 - u[0] * u[0]) * u[1] - u[0]);
187c4762a1bSJed Brown 
1889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
191c4762a1bSJed Brown   PetscFunctionReturn(0);
192c4762a1bSJed Brown }
193c4762a1bSJed Brown 
1949371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) {
195c4762a1bSJed Brown   User               user     = (User)ctx;
196c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
197c4762a1bSJed Brown   PetscScalar        J[2][2];
198c4762a1bSJed Brown   const PetscScalar *u;
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   PetscFunctionBeginUser;
2019566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
202c4762a1bSJed Brown 
2039371c9d4SSatish Balay   J[0][0] = a;
2049371c9d4SSatish Balay   J[0][1] = -1.0;
2059371c9d4SSatish Balay   J[1][0] = user->mu * (1.0 + 2.0 * u[0] * u[1]);
2069371c9d4SSatish Balay   J[1][1] = a - user->mu * (1.0 - u[0] * u[0]);
2079566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
2089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2099566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2109566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
211c4762a1bSJed Brown   if (A != B) {
2129566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2139566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
214c4762a1bSJed Brown   }
215c4762a1bSJed Brown 
2169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
217c4762a1bSJed Brown   PetscFunctionReturn(0);
218c4762a1bSJed Brown }
219c4762a1bSJed Brown 
2209371c9d4SSatish Balay static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec U, Mat A, void *ctx) {
221c4762a1bSJed Brown   PetscInt           row[] = {0, 1}, col[] = {0};
222c4762a1bSJed Brown   PetscScalar        J[2][1];
223c4762a1bSJed Brown   const PetscScalar *u;
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   PetscFunctionBeginUser;
2269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   J[0][0] = 0;
229c4762a1bSJed Brown   J[1][0] = (1. - u[0] * u[0]) * u[1] - u[0];
2309566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
2319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
234c4762a1bSJed Brown 
2359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
236c4762a1bSJed Brown   PetscFunctionReturn(0);
237c4762a1bSJed Brown }
238c4762a1bSJed Brown 
2399371c9d4SSatish Balay static PetscErrorCode IHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
240c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
241c4762a1bSJed Brown   PetscScalar       *vhv;
242c4762a1bSJed Brown   PetscScalar        dJdU[2][2][2] = {{{0}}};
243c4762a1bSJed Brown   PetscInt           i, j, k;
244c4762a1bSJed Brown   User               user = (User)ctx;
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   PetscFunctionBeginUser;
2479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
2489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
2499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
2509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   dJdU[1][0][0] = 2. * user->mu * u[1];
253c4762a1bSJed Brown   dJdU[1][1][0] = 2. * user->mu * u[0];
254c4762a1bSJed Brown   dJdU[1][0][1] = 2. * user->mu * u[0];
255c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
256c4762a1bSJed Brown     vhv[j] = 0;
257c4762a1bSJed Brown     for (k = 0; k < 2; k++)
2589371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
259c4762a1bSJed Brown   }
260c4762a1bSJed Brown 
2619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
2639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
2649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
265c4762a1bSJed Brown   PetscFunctionReturn(0);
266c4762a1bSJed Brown }
267c4762a1bSJed Brown 
2689371c9d4SSatish Balay static PetscErrorCode IHessianProductUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
269c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
270c4762a1bSJed Brown   PetscScalar       *vhv;
271c4762a1bSJed Brown   PetscScalar        dJdP[2][2][1] = {{{0}}};
272c4762a1bSJed Brown   PetscInt           i, j, k;
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   PetscFunctionBeginUser;
2759566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
2769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
2779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
2789566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   dJdP[1][0][0] = 1. + 2. * u[0] * u[1];
281c4762a1bSJed Brown   dJdP[1][1][0] = u[0] * u[0] - 1.;
282c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
283c4762a1bSJed Brown     vhv[j] = 0;
284c4762a1bSJed Brown     for (k = 0; k < 1; k++)
2859371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdP[i][j][k] * vr[k];
286c4762a1bSJed Brown   }
287c4762a1bSJed Brown 
2889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
2909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
2919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
292c4762a1bSJed Brown   PetscFunctionReturn(0);
293c4762a1bSJed Brown }
294c4762a1bSJed Brown 
2959371c9d4SSatish Balay static PetscErrorCode IHessianProductPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
296c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
297c4762a1bSJed Brown   PetscScalar       *vhv;
298c4762a1bSJed Brown   PetscScalar        dJdU[2][1][2] = {{{0}}};
299c4762a1bSJed Brown   PetscInt           i, j, k;
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   PetscFunctionBeginUser;
3029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
3039566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
3049566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
3059566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   dJdU[1][0][0] = 1. + 2. * u[1] * u[0];
308c4762a1bSJed Brown   dJdU[1][0][1] = u[0] * u[0] - 1.;
309c4762a1bSJed Brown   for (j = 0; j < 1; j++) {
310c4762a1bSJed Brown     vhv[j] = 0;
311c4762a1bSJed Brown     for (k = 0; k < 2; k++)
3129371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
313c4762a1bSJed Brown   }
314c4762a1bSJed Brown 
3159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
3169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
3179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
3189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
319c4762a1bSJed Brown   PetscFunctionReturn(0);
320c4762a1bSJed Brown }
321c4762a1bSJed Brown 
3229371c9d4SSatish Balay static PetscErrorCode IHessianProductPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
323c4762a1bSJed Brown   PetscFunctionBeginUser;
324c4762a1bSJed Brown   PetscFunctionReturn(0);
325c4762a1bSJed Brown }
326c4762a1bSJed Brown 
327c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
3289371c9d4SSatish Balay static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) {
329c4762a1bSJed Brown   const PetscScalar *x;
330c4762a1bSJed Brown   PetscReal          tfinal, dt;
331c4762a1bSJed Brown   User               user = (User)ctx;
332c4762a1bSJed Brown   Vec                interpolatedX;
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   PetscFunctionBeginUser;
3359566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
3369566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts, &tfinal));
337c4762a1bSJed Brown 
338c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
3399566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &interpolatedX));
3409566063dSJacob Faibussowitsch     PetscCall(TSInterpolate(ts, user->next_output, interpolatedX));
3419566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(interpolatedX, &x));
3429371c9d4SSatish 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(x[0]), (double)PetscRealPart(x[1])));
3439566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(interpolatedX, &x));
3449566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&interpolatedX));
345c4762a1bSJed Brown     user->next_output += PetscRealConstant(0.1);
346c4762a1bSJed Brown   }
347c4762a1bSJed Brown   PetscFunctionReturn(0);
348c4762a1bSJed Brown }
349c4762a1bSJed Brown 
3509371c9d4SSatish Balay int main(int argc, char **argv) {
351c4762a1bSJed Brown   Vec                P;
352c4762a1bSJed Brown   PetscBool          monitor = PETSC_FALSE;
353c4762a1bSJed Brown   PetscScalar       *x_ptr;
354c4762a1bSJed Brown   const PetscScalar *y_ptr;
355c4762a1bSJed Brown   PetscMPIInt        size;
356c4762a1bSJed Brown   struct _n_User     user;
357c4762a1bSJed Brown   Tao                tao;
358c4762a1bSJed Brown   KSP                ksp;
359c4762a1bSJed Brown   PC                 pc;
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
362c4762a1bSJed Brown      Initialize program
363c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
364327415f7SBarry Smith   PetscFunctionBeginUser;
3659566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
3673c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
368c4762a1bSJed Brown 
369c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
3709566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
3719566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOBQNLS));
372c4762a1bSJed Brown 
373c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
374c4762a1bSJed Brown     Set runtime options
375c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
376c4762a1bSJed Brown   user.next_output  = 0.0;
377c4762a1bSJed Brown   user.mu           = PetscRealConstant(1.0e3);
378c4762a1bSJed Brown   user.ftime        = PetscRealConstant(0.5);
379c4762a1bSJed Brown   user.implicitform = PETSC_TRUE;
380c4762a1bSJed Brown 
3819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
3829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
3839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &user.implicitform, NULL));
384c4762a1bSJed Brown 
385c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
3869566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.A));
3879566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
3889566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.A));
3899566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.A));
3909566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.U, NULL));
3919566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Lambda[0], NULL));
3929566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Lambda2[0], NULL));
3939566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Ihp1[0], NULL));
3949566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Ihp2[0], NULL));
395c4762a1bSJed Brown 
3969566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Jacp));
3979566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
3989566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.Jacp));
3999566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.Jacp));
4009566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &user.Dir, NULL));
4019566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &user.Mup[0], NULL));
4029566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &user.Mup2[0], NULL));
4039566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &user.Ihp3[0], NULL));
4049566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &user.Ihp4[0], NULL));
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   /* Create timestepping solver context */
4079566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &user.ts));
4089566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(user.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
409c4762a1bSJed Brown   if (user.implicitform) {
4109566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(user.ts, NULL, IFunction, &user));
4119566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(user.ts, user.A, user.A, IJacobian, &user));
4129566063dSJacob Faibussowitsch     PetscCall(TSSetType(user.ts, TSCN));
413c4762a1bSJed Brown   } else {
4149566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(user.ts, NULL, RHSFunction, &user));
4159566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(user.ts, user.A, user.A, RHSJacobian, &user));
4169566063dSJacob Faibussowitsch     PetscCall(TSSetType(user.ts, TSRK));
417c4762a1bSJed Brown   }
4189566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(user.ts, user.Jacp, RHSJacobianP, &user));
4199566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(user.ts, user.ftime));
4209566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(user.ts, TS_EXACTFINALTIME_MATCHSTEP));
421c4762a1bSJed Brown 
422*48a46eb9SPierre Jolivet   if (monitor) PetscCall(TSMonitorSet(user.ts, Monitor, &user, NULL));
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   /* Set ODE initial conditions */
4259566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.U, &x_ptr));
426c4762a1bSJed Brown   x_ptr[0] = 2.0;
427c4762a1bSJed Brown   x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu);
4289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.U, &x_ptr));
4299566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(user.ts, PetscRealConstant(0.001)));
430c4762a1bSJed Brown 
431c4762a1bSJed Brown   /* Set runtime options */
4329566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(user.ts));
433c4762a1bSJed Brown 
4349566063dSJacob Faibussowitsch   PetscCall(TSSolve(user.ts, user.U));
4359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user.U, &y_ptr));
436c4762a1bSJed Brown   user.ob[0] = y_ptr[0];
437c4762a1bSJed Brown   user.ob[1] = y_ptr[1];
4389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user.U, &y_ptr));
439c4762a1bSJed Brown 
440c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used.
441c4762a1bSJed Brown      Skip checkpointing for the first TSSolve since no adjoint run follows it.
442c4762a1bSJed Brown    */
4439566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(user.ts));
444c4762a1bSJed Brown 
445c4762a1bSJed Brown   /* Optimization starts */
4469566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.H));
4479566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.H, PETSC_DECIDE, PETSC_DECIDE, 1, 1));
4489566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.H)); /* Hessian should be symmetric. Do we need to do MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE) ? */
449c4762a1bSJed Brown 
450c4762a1bSJed Brown   /* Set initial solution guess */
4519566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &P, NULL));
4529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(P, &x_ptr));
453c4762a1bSJed Brown   x_ptr[0] = PetscRealConstant(1.2);
4549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(P, &x_ptr));
4559566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, P));
456c4762a1bSJed Brown 
457c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
4589566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
4599566063dSJacob Faibussowitsch   PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user));
460c4762a1bSJed Brown 
461c4762a1bSJed Brown   /* Check for any TAO command line options */
4629566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao, &ksp));
463c4762a1bSJed Brown   if (ksp) {
4649566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
4659566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, PCNONE));
466c4762a1bSJed Brown   }
4679566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
468c4762a1bSJed Brown 
4699566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
470c4762a1bSJed Brown 
4719566063dSJacob Faibussowitsch   PetscCall(VecView(P, PETSC_VIEWER_STDOUT_WORLD));
472c4762a1bSJed Brown   /* Free TAO data structures */
4739566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
474c4762a1bSJed Brown 
475c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
476c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
477c4762a1bSJed Brown      are no longer needed.
478c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4799566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.H));
4809566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.A));
4819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.U));
4829566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Jacp));
4839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Lambda[0]));
4849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Mup[0]));
4859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Lambda2[0]));
4869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Mup2[0]));
4879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp1[0]));
4889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp2[0]));
4899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp3[0]));
4909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp4[0]));
4919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Dir));
4929566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&user.ts));
4939566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&P));
4949566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
495b122ec5aSJacob Faibussowitsch   return 0;
496c4762a1bSJed Brown }
497c4762a1bSJed Brown 
498c4762a1bSJed Brown /* ------------------------------------------------------------------ */
499c4762a1bSJed Brown /*
500c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
501c4762a1bSJed Brown 
502c4762a1bSJed Brown    Input Parameters:
503c4762a1bSJed Brown    tao - the Tao context
504c4762a1bSJed Brown    X   - the input vector
505a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
506c4762a1bSJed Brown 
507c4762a1bSJed Brown    Output Parameters:
508c4762a1bSJed Brown    f   - the newly evaluated function
509c4762a1bSJed Brown    G   - the newly evaluated gradient
510c4762a1bSJed Brown */
5119371c9d4SSatish Balay PetscErrorCode FormFunctionGradient(Tao tao, Vec P, PetscReal *f, Vec G, void *ctx) {
512c4762a1bSJed Brown   User               user_ptr = (User)ctx;
513c4762a1bSJed Brown   TS                 ts       = user_ptr->ts;
514c4762a1bSJed Brown   PetscScalar       *x_ptr, *g;
515c4762a1bSJed Brown   const PetscScalar *y_ptr;
516c4762a1bSJed Brown 
517c4762a1bSJed Brown   PetscFunctionBeginUser;
5189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P, &y_ptr));
519c4762a1bSJed Brown   user_ptr->mu = y_ptr[0];
5209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P, &y_ptr));
521c4762a1bSJed Brown 
5229566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
5239566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
5249566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, PetscRealConstant(0.001))); /* can be overwritten by command line options */
5259566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
5269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->U, &x_ptr));
527c4762a1bSJed Brown   x_ptr[0] = 2.0;
528c4762a1bSJed Brown   x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user_ptr->mu) - 292.0 / (2187.0 * user_ptr->mu * user_ptr->mu);
5299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->U, &x_ptr));
530c4762a1bSJed Brown 
5319566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, user_ptr->U));
532c4762a1bSJed Brown 
5339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user_ptr->U, &y_ptr));
534c4762a1bSJed Brown   *f = (y_ptr[0] - user_ptr->ob[0]) * (y_ptr[0] - user_ptr->ob[0]) + (y_ptr[1] - user_ptr->ob[1]) * (y_ptr[1] - user_ptr->ob[1]);
535c4762a1bSJed Brown 
536c4762a1bSJed Brown   /*   Reset initial conditions for the adjoint integration */
5379566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Lambda[0], &x_ptr));
538c4762a1bSJed Brown   x_ptr[0] = 2. * (y_ptr[0] - user_ptr->ob[0]);
539c4762a1bSJed Brown   x_ptr[1] = 2. * (y_ptr[1] - user_ptr->ob[1]);
5409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user_ptr->U, &y_ptr));
5419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Lambda[0], &x_ptr));
542c4762a1bSJed Brown 
5439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Mup[0], &x_ptr));
544c4762a1bSJed Brown   x_ptr[0] = 0.0;
5459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Mup[0], &x_ptr));
5469566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, user_ptr->Lambda, user_ptr->Mup));
547c4762a1bSJed Brown 
5489566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
549c4762a1bSJed Brown 
5509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Mup[0], &x_ptr));
5519566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user_ptr->Lambda[0], &y_ptr));
5529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G, &g));
553c4762a1bSJed Brown   g[0] = y_ptr[1] * (-10.0 / (81.0 * user_ptr->mu * user_ptr->mu) + 2.0 * 292.0 / (2187.0 * user_ptr->mu * user_ptr->mu * user_ptr->mu)) + x_ptr[0];
5549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Mup[0], &x_ptr));
5559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user_ptr->Lambda[0], &y_ptr));
5569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G, &g));
557c4762a1bSJed Brown   PetscFunctionReturn(0);
558c4762a1bSJed Brown }
559c4762a1bSJed Brown 
5609371c9d4SSatish Balay PetscErrorCode FormHessian(Tao tao, Vec P, Mat H, Mat Hpre, void *ctx) {
561c4762a1bSJed Brown   User           user_ptr = (User)ctx;
562c4762a1bSJed Brown   PetscScalar    harr[1];
563c4762a1bSJed Brown   const PetscInt rows[1] = {0};
564c4762a1bSJed Brown   PetscInt       col     = 0;
565c4762a1bSJed Brown 
566c4762a1bSJed Brown   PetscFunctionBeginUser;
5679566063dSJacob Faibussowitsch   PetscCall(Adjoint2(P, harr, user_ptr));
5689566063dSJacob Faibussowitsch   PetscCall(MatSetValues(H, 1, rows, 1, &col, harr, INSERT_VALUES));
569c4762a1bSJed Brown 
5709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
5719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
572c4762a1bSJed Brown   if (H != Hpre) {
5739566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Hpre, MAT_FINAL_ASSEMBLY));
5749566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Hpre, MAT_FINAL_ASSEMBLY));
575c4762a1bSJed Brown   }
576c4762a1bSJed Brown   PetscFunctionReturn(0);
577c4762a1bSJed Brown }
578c4762a1bSJed Brown 
5799371c9d4SSatish Balay PetscErrorCode Adjoint2(Vec P, PetscScalar arr[], User ctx) {
580c4762a1bSJed Brown   TS                 ts = ctx->ts;
581c4762a1bSJed Brown   const PetscScalar *z_ptr;
582c4762a1bSJed Brown   PetscScalar       *x_ptr, *y_ptr, dzdp, dzdp2;
583c4762a1bSJed Brown   Mat                tlmsen;
584c4762a1bSJed Brown 
585c4762a1bSJed Brown   PetscFunctionBeginUser;
586c4762a1bSJed Brown   /* Reset TSAdjoint so that AdjointSetUp will be called again */
5879566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
588c4762a1bSJed Brown 
589c4762a1bSJed Brown   /* The directional vector should be 1 since it is one-dimensional */
5909566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Dir, &x_ptr));
591c4762a1bSJed Brown   x_ptr[0] = 1.;
5929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Dir, &x_ptr));
593c4762a1bSJed Brown 
5949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P, &z_ptr));
595c4762a1bSJed Brown   ctx->mu = z_ptr[0];
5969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P, &z_ptr));
597c4762a1bSJed Brown 
598c4762a1bSJed Brown   dzdp  = -10.0 / (81.0 * ctx->mu * ctx->mu) + 2.0 * 292.0 / (2187.0 * ctx->mu * ctx->mu * ctx->mu);
599c4762a1bSJed Brown   dzdp2 = 2. * 10.0 / (81.0 * ctx->mu * ctx->mu * ctx->mu) - 3.0 * 2.0 * 292.0 / (2187.0 * ctx->mu * ctx->mu * ctx->mu * ctx->mu);
600c4762a1bSJed Brown 
6019566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
6029566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
6039566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, PetscRealConstant(0.001))); /* can be overwritten by command line options */
6049566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
6059566063dSJacob Faibussowitsch   PetscCall(TSSetCostHessianProducts(ts, 1, ctx->Lambda2, ctx->Mup2, ctx->Dir));
606c4762a1bSJed Brown 
6079566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(ctx->Jacp));
6089566063dSJacob Faibussowitsch   PetscCall(MatSetValue(ctx->Jacp, 1, 0, dzdp, INSERT_VALUES));
6099566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(ctx->Jacp, MAT_FINAL_ASSEMBLY));
6109566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(ctx->Jacp, MAT_FINAL_ASSEMBLY));
611c4762a1bSJed Brown 
6129566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetForward(ts, ctx->Jacp));
6139566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->U, &y_ptr));
614c4762a1bSJed Brown   y_ptr[0] = 2.0;
615c4762a1bSJed Brown   y_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * ctx->mu) - 292.0 / (2187.0 * ctx->mu * ctx->mu);
6169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->U, &y_ptr));
6179566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, ctx->U));
618c4762a1bSJed Brown 
619c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
6209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->U, &z_ptr));
6219566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda[0], &y_ptr));
622c4762a1bSJed Brown   y_ptr[0] = 2. * (z_ptr[0] - ctx->ob[0]);
623c4762a1bSJed Brown   y_ptr[1] = 2. * (z_ptr[1] - ctx->ob[1]);
6249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda[0], &y_ptr));
6259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->U, &z_ptr));
6269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Mup[0], &y_ptr));
627c4762a1bSJed Brown   y_ptr[0] = 0.0;
6289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Mup[0], &y_ptr));
6299566063dSJacob Faibussowitsch   PetscCall(TSForwardGetSensitivities(ts, NULL, &tlmsen));
6309566063dSJacob Faibussowitsch   PetscCall(MatDenseGetColumn(tlmsen, 0, &x_ptr));
6319566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda2[0], &y_ptr));
632c4762a1bSJed Brown   y_ptr[0] = 2. * x_ptr[0];
633c4762a1bSJed Brown   y_ptr[1] = 2. * x_ptr[1];
6349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0], &y_ptr));
6359566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Mup2[0], &y_ptr));
636c4762a1bSJed Brown   y_ptr[0] = 0.0;
6379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Mup2[0], &y_ptr));
6389566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreColumn(tlmsen, &x_ptr));
6399566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, ctx->Lambda, ctx->Mup));
640c4762a1bSJed Brown   if (ctx->implicitform) {
6419566063dSJacob Faibussowitsch     PetscCall(TSSetIHessianProduct(ts, ctx->Ihp1, IHessianProductUU, ctx->Ihp2, IHessianProductUP, ctx->Ihp3, IHessianProductPU, ctx->Ihp4, IHessianProductPP, ctx));
642c4762a1bSJed Brown   } else {
6439566063dSJacob Faibussowitsch     PetscCall(TSSetRHSHessianProduct(ts, ctx->Ihp1, RHSHessianProductUU, ctx->Ihp2, RHSHessianProductUP, ctx->Ihp3, RHSHessianProductPU, ctx->Ihp4, RHSHessianProductPP, ctx));
644c4762a1bSJed Brown   }
6459566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
646c4762a1bSJed Brown 
6479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda[0], &x_ptr));
6489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda2[0], &y_ptr));
6499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->Mup2[0], &z_ptr));
650c4762a1bSJed Brown 
651c4762a1bSJed Brown   arr[0] = x_ptr[1] * dzdp2 + y_ptr[1] * dzdp2 + z_ptr[0];
652c4762a1bSJed Brown 
6539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0], &x_ptr));
6549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0], &y_ptr));
6559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->Mup2[0], &z_ptr));
656c4762a1bSJed Brown 
657c4762a1bSJed Brown   /* Disable second-order adjoint mode */
6589566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
6599566063dSJacob Faibussowitsch   PetscCall(TSAdjointResetForward(ts));
660c4762a1bSJed Brown   PetscFunctionReturn(0);
661c4762a1bSJed Brown }
662c4762a1bSJed Brown 
663c4762a1bSJed Brown /*TEST
664c4762a1bSJed Brown     build:
665c4762a1bSJed Brown       requires: !complex !single
666c4762a1bSJed Brown     test:
667c4762a1bSJed Brown       args:  -implicitform 0 -ts_type rk -ts_adapt_type none -mu 10 -ts_dt 0.1 -viewer_binary_skip_info -tao_monitor -tao_view
668c4762a1bSJed Brown       output_file: output/ex20opt_p_1.out
669c4762a1bSJed Brown 
670c4762a1bSJed Brown     test:
671c4762a1bSJed Brown       suffix: 2
672c4762a1bSJed Brown       args:  -implicitform 0 -ts_type rk -ts_adapt_type none -mu 10 -ts_dt 0.01 -viewer_binary_skip_info -tao_monitor -tao_type bntr -tao_bnk_pc_type none
673c4762a1bSJed Brown       output_file: output/ex20opt_p_2.out
674c4762a1bSJed Brown 
675c4762a1bSJed Brown     test:
676c4762a1bSJed Brown       suffix: 3
677c4762a1bSJed Brown       args:  -ts_type cn -ts_adapt_type none -mu 100 -ts_dt 0.01 -viewer_binary_skip_info -tao_monitor -tao_view
678c4762a1bSJed Brown       output_file: output/ex20opt_p_3.out
679c4762a1bSJed Brown 
680c4762a1bSJed Brown     test:
681c4762a1bSJed Brown       suffix: 4
682c4762a1bSJed Brown       args:  -ts_type cn -ts_adapt_type none -mu 100 -ts_dt 0.01 -viewer_binary_skip_info -tao_monitor -tao_type bntr -tao_bnk_pc_type none
683c4762a1bSJed Brown       output_file: output/ex20opt_p_4.out
684c4762a1bSJed Brown 
685c4762a1bSJed Brown TEST*/
686