xref: /petsc/src/ts/tutorials/ex20opt_p.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation.\n\
2c4762a1bSJed Brown Input parameters include:\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /* ------------------------------------------------------------------------
5c4762a1bSJed Brown 
6c4762a1bSJed Brown   Notes:
7c4762a1bSJed Brown   This code demonstrates how to solve a DAE-constrained optimization problem with TAO, TSAdjoint and TS.
8c4762a1bSJed Brown   The nonlinear problem is written in a DAE equivalent form.
9c4762a1bSJed Brown   The objective is to minimize the difference between observation and model prediction by finding an optimal value for parameter \mu.
10c4762a1bSJed Brown   The gradient is computed with the discrete adjoint of an implicit theta method, see ex20adj.c for details.
11c4762a1bSJed Brown   ------------------------------------------------------------------------- */
12c4762a1bSJed Brown #include <petsctao.h>
13c4762a1bSJed Brown #include <petscts.h>
14c4762a1bSJed Brown 
15c4762a1bSJed Brown typedef struct _n_User *User;
16c4762a1bSJed Brown struct _n_User {
17c4762a1bSJed Brown   TS        ts;
18c4762a1bSJed Brown   PetscReal mu;
19c4762a1bSJed Brown   PetscReal next_output;
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   /* Sensitivity analysis support */
22c4762a1bSJed Brown   PetscReal ftime;
23c4762a1bSJed Brown   Mat       A;                    /* Jacobian matrix */
24c4762a1bSJed Brown   Mat       Jacp;                 /* JacobianP matrix */
25c4762a1bSJed Brown   Mat       H;                    /* Hessian matrix for optimization */
26c4762a1bSJed Brown   Vec       U, Lambda[1], Mup[1]; /* adjoint variables */
27c4762a1bSJed Brown   Vec       Lambda2[1], Mup2[1];  /* second-order adjoint variables */
28c4762a1bSJed Brown   Vec       Ihp1[1];              /* working space for Hessian evaluations */
29c4762a1bSJed Brown   Vec       Ihp2[1];              /* working space for Hessian evaluations */
30c4762a1bSJed Brown   Vec       Ihp3[1];              /* working space for Hessian evaluations */
31c4762a1bSJed Brown   Vec       Ihp4[1];              /* working space for Hessian evaluations */
32c4762a1bSJed Brown   Vec       Dir;                  /* direction vector */
33c4762a1bSJed Brown   PetscReal ob[2];                /* observation used by the cost function */
34c4762a1bSJed Brown   PetscBool implicitform;         /* implicit ODE? */
35c4762a1bSJed Brown };
36c4762a1bSJed Brown 
37c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
38c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
39c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec, PetscScalar[], User);
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE  -------------------- */
42c4762a1bSJed Brown 
RHSFunction(TS ts,PetscReal t,Vec U,Vec F,PetscCtx ctx)43*2a8381b2SBarry Smith static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, PetscCtx ctx)
44d71ae5a4SJacob Faibussowitsch {
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));
563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,PetscCtx ctx)59*2a8381b2SBarry Smith static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, PetscCtx ctx)
60d71ae5a4SJacob Faibussowitsch {
61c4762a1bSJed Brown   User               user     = (User)ctx;
62c4762a1bSJed Brown   PetscReal          mu       = user->mu;
63c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
64c4762a1bSJed Brown   PetscScalar        J[2][2];
65c4762a1bSJed Brown   const PetscScalar *u;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   PetscFunctionBeginUser;
689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   J[0][0] = 0;
71c4762a1bSJed Brown   J[1][0] = -mu * (2.0 * u[1] * u[0] + 1.);
72c4762a1bSJed Brown   J[0][1] = 1.0;
73c4762a1bSJed Brown   J[1][1] = mu * (1.0 - u[0] * u[0]);
749566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
77c4762a1bSJed Brown   if (B && A != B) {
789566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
799566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
80c4762a1bSJed Brown   }
81c4762a1bSJed Brown 
829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
84c4762a1bSJed Brown }
85c4762a1bSJed Brown 
RHSHessianProductUU(TS ts,PetscReal t,Vec U,Vec * Vl,Vec Vr,Vec * VHV,PetscCtx ctx)86*2a8381b2SBarry Smith static PetscErrorCode RHSHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, PetscCtx ctx)
87d71ae5a4SJacob Faibussowitsch {
88c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
89c4762a1bSJed Brown   PetscScalar       *vhv;
90c4762a1bSJed Brown   PetscScalar        dJdU[2][2][2] = {{{0}}};
91c4762a1bSJed Brown   PetscInt           i, j, k;
92c4762a1bSJed Brown   User               user = (User)ctx;
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   PetscFunctionBeginUser;
959566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
969566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
979566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
989566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   dJdU[1][0][0] = -2. * user->mu * u[1];
101c4762a1bSJed Brown   dJdU[1][1][0] = -2. * user->mu * u[0];
102c4762a1bSJed Brown   dJdU[1][0][1] = -2. * user->mu * u[0];
103c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
104c4762a1bSJed Brown     vhv[j] = 0;
105c4762a1bSJed Brown     for (k = 0; k < 2; k++)
1069371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
107c4762a1bSJed Brown   }
108c4762a1bSJed Brown 
1099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
1119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
RHSHessianProductUP(TS ts,PetscReal t,Vec U,Vec * Vl,Vec Vr,Vec * VHV,PetscCtx ctx)116*2a8381b2SBarry Smith static PetscErrorCode RHSHessianProductUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, PetscCtx ctx)
117d71ae5a4SJacob Faibussowitsch {
118c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
119c4762a1bSJed Brown   PetscScalar       *vhv;
120c4762a1bSJed Brown   PetscScalar        dJdP[2][2][1] = {{{0}}};
121c4762a1bSJed Brown   PetscInt           i, j, k;
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   PetscFunctionBeginUser;
1249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
1269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
1279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   dJdP[1][0][0] = -(1. + 2. * u[0] * u[1]);
130c4762a1bSJed Brown   dJdP[1][1][0] = 1. - u[0] * u[0];
131c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
132c4762a1bSJed Brown     vhv[j] = 0;
133c4762a1bSJed Brown     for (k = 0; k < 1; k++)
1349371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdP[i][j][k] * vr[k];
135c4762a1bSJed Brown   }
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
1399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
1413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
142c4762a1bSJed Brown }
143c4762a1bSJed Brown 
RHSHessianProductPU(TS ts,PetscReal t,Vec U,Vec * Vl,Vec Vr,Vec * VHV,PetscCtx ctx)144*2a8381b2SBarry Smith static PetscErrorCode RHSHessianProductPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, PetscCtx ctx)
145d71ae5a4SJacob Faibussowitsch {
146c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
147c4762a1bSJed Brown   PetscScalar       *vhv;
148c4762a1bSJed Brown   PetscScalar        dJdU[2][1][2] = {{{0}}};
149c4762a1bSJed Brown   PetscInt           i, j, k;
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   PetscFunctionBeginUser;
1529566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
1549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
1559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   dJdU[1][0][0] = -1. - 2. * u[1] * u[0];
158c4762a1bSJed Brown   dJdU[1][0][1] = 1. - u[0] * u[0];
159c4762a1bSJed Brown   for (j = 0; j < 1; j++) {
160c4762a1bSJed Brown     vhv[j] = 0;
161c4762a1bSJed Brown     for (k = 0; k < 2; k++)
1629371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
163c4762a1bSJed Brown   }
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
1679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170c4762a1bSJed Brown }
171c4762a1bSJed Brown 
RHSHessianProductPP(TS ts,PetscReal t,Vec U,Vec * Vl,Vec Vr,Vec * VHV,PetscCtx ctx)172*2a8381b2SBarry Smith static PetscErrorCode RHSHessianProductPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, PetscCtx ctx)
173d71ae5a4SJacob Faibussowitsch {
174c4762a1bSJed Brown   PetscFunctionBeginUser;
1753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
176c4762a1bSJed Brown }
177c4762a1bSJed Brown 
178c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE  -------------------- */
179c4762a1bSJed Brown 
IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,PetscCtx ctx)180*2a8381b2SBarry Smith static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx)
181d71ae5a4SJacob Faibussowitsch {
182c4762a1bSJed Brown   User               user = (User)ctx;
183c4762a1bSJed Brown   PetscScalar       *f;
184c4762a1bSJed Brown   const PetscScalar *u, *udot;
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   PetscFunctionBeginUser;
1879566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   f[0] = udot[0] - u[1];
192c4762a1bSJed Brown   f[1] = udot[1] - user->mu * ((1.0 - u[0] * u[0]) * u[1] - u[0]);
193c4762a1bSJed Brown 
1949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
198c4762a1bSJed Brown }
199c4762a1bSJed Brown 
IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,PetscCtx ctx)200*2a8381b2SBarry Smith static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
201d71ae5a4SJacob Faibussowitsch {
202c4762a1bSJed Brown   User               user     = (User)ctx;
203c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
204c4762a1bSJed Brown   PetscScalar        J[2][2];
205c4762a1bSJed Brown   const PetscScalar *u;
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   PetscFunctionBeginUser;
2089566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
209c4762a1bSJed Brown 
2109371c9d4SSatish Balay   J[0][0] = a;
2119371c9d4SSatish Balay   J[0][1] = -1.0;
2129371c9d4SSatish Balay   J[1][0] = user->mu * (1.0 + 2.0 * u[0] * u[1]);
2139371c9d4SSatish Balay   J[1][1] = a - user->mu * (1.0 - u[0] * u[0]);
2149566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
2159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2169566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2179566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
218c4762a1bSJed Brown   if (A != B) {
2199566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2209566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
221c4762a1bSJed Brown   }
222c4762a1bSJed Brown 
2239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
225c4762a1bSJed Brown }
226c4762a1bSJed Brown 
RHSJacobianP(TS ts,PetscReal t,Vec U,Mat A,PetscCtx ctx)227*2a8381b2SBarry Smith static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec U, Mat A, PetscCtx ctx)
228d71ae5a4SJacob Faibussowitsch {
229c4762a1bSJed Brown   PetscInt           row[] = {0, 1}, col[] = {0};
230c4762a1bSJed Brown   PetscScalar        J[2][1];
231c4762a1bSJed Brown   const PetscScalar *u;
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   PetscFunctionBeginUser;
2349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   J[0][0] = 0;
237c4762a1bSJed Brown   J[1][0] = (1. - u[0] * u[0]) * u[1] - u[0];
2389566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
2399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2419566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
242c4762a1bSJed Brown 
2439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
245c4762a1bSJed Brown }
246c4762a1bSJed Brown 
IHessianProductUU(TS ts,PetscReal t,Vec U,Vec * Vl,Vec Vr,Vec * VHV,PetscCtx ctx)247*2a8381b2SBarry Smith static PetscErrorCode IHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, PetscCtx ctx)
248d71ae5a4SJacob Faibussowitsch {
249c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
250c4762a1bSJed Brown   PetscScalar       *vhv;
251c4762a1bSJed Brown   PetscScalar        dJdU[2][2][2] = {{{0}}};
252c4762a1bSJed Brown   PetscInt           i, j, k;
253c4762a1bSJed Brown   User               user = (User)ctx;
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   PetscFunctionBeginUser;
2569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
2579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
2589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
2599566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   dJdU[1][0][0] = 2. * user->mu * u[1];
262c4762a1bSJed Brown   dJdU[1][1][0] = 2. * user->mu * u[0];
263c4762a1bSJed Brown   dJdU[1][0][1] = 2. * user->mu * u[0];
264c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
265c4762a1bSJed Brown     vhv[j] = 0;
266c4762a1bSJed Brown     for (k = 0; k < 2; k++)
2679371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
268c4762a1bSJed Brown   }
269c4762a1bSJed Brown 
2709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
2729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
2739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275c4762a1bSJed Brown }
276c4762a1bSJed Brown 
IHessianProductUP(TS ts,PetscReal t,Vec U,Vec * Vl,Vec Vr,Vec * VHV,PetscCtx ctx)277*2a8381b2SBarry Smith static PetscErrorCode IHessianProductUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, PetscCtx ctx)
278d71ae5a4SJacob Faibussowitsch {
279c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
280c4762a1bSJed Brown   PetscScalar       *vhv;
281c4762a1bSJed Brown   PetscScalar        dJdP[2][2][1] = {{{0}}};
282c4762a1bSJed Brown   PetscInt           i, j, k;
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   PetscFunctionBeginUser;
2859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
2869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
2879566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
2889566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   dJdP[1][0][0] = 1. + 2. * u[0] * u[1];
291c4762a1bSJed Brown   dJdP[1][1][0] = u[0] * u[0] - 1.;
292c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
293c4762a1bSJed Brown     vhv[j] = 0;
294c4762a1bSJed Brown     for (k = 0; k < 1; k++)
2959371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdP[i][j][k] * vr[k];
296c4762a1bSJed Brown   }
297c4762a1bSJed Brown 
2989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
3009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
3019566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
3023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
303c4762a1bSJed Brown }
304c4762a1bSJed Brown 
IHessianProductPU(TS ts,PetscReal t,Vec U,Vec * Vl,Vec Vr,Vec * VHV,PetscCtx ctx)305*2a8381b2SBarry Smith static PetscErrorCode IHessianProductPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, PetscCtx ctx)
306d71ae5a4SJacob Faibussowitsch {
307c4762a1bSJed Brown   const PetscScalar *vl, *vr, *u;
308c4762a1bSJed Brown   PetscScalar       *vhv;
309c4762a1bSJed Brown   PetscScalar        dJdU[2][1][2] = {{{0}}};
310c4762a1bSJed Brown   PetscInt           i, j, k;
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   PetscFunctionBeginUser;
3139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
3149566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
3159566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
3169566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   dJdU[1][0][0] = 1. + 2. * u[1] * u[0];
319c4762a1bSJed Brown   dJdU[1][0][1] = u[0] * u[0] - 1.;
320c4762a1bSJed Brown   for (j = 0; j < 1; j++) {
321c4762a1bSJed Brown     vhv[j] = 0;
322c4762a1bSJed Brown     for (k = 0; k < 2; k++)
3239371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k];
324c4762a1bSJed Brown   }
325c4762a1bSJed Brown 
3269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
3279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
3289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
3299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
3303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
331c4762a1bSJed Brown }
332c4762a1bSJed Brown 
IHessianProductPP(TS ts,PetscReal t,Vec U,Vec * Vl,Vec Vr,Vec * VHV,PetscCtx ctx)333*2a8381b2SBarry Smith static PetscErrorCode IHessianProductPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, PetscCtx ctx)
334d71ae5a4SJacob Faibussowitsch {
335c4762a1bSJed Brown   PetscFunctionBeginUser;
3363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
337c4762a1bSJed Brown }
338c4762a1bSJed Brown 
339c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
Monitor(TS ts,PetscInt step,PetscReal t,Vec X,PetscCtx ctx)340*2a8381b2SBarry Smith static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, PetscCtx ctx)
341d71ae5a4SJacob Faibussowitsch {
342c4762a1bSJed Brown   const PetscScalar *x;
343c4762a1bSJed Brown   PetscReal          tfinal, dt;
344c4762a1bSJed Brown   User               user = (User)ctx;
345c4762a1bSJed Brown   Vec                interpolatedX;
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   PetscFunctionBeginUser;
3489566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
3499566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts, &tfinal));
350c4762a1bSJed Brown 
351c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
3529566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &interpolatedX));
3539566063dSJacob Faibussowitsch     PetscCall(TSInterpolate(ts, user->next_output, interpolatedX));
3549566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(interpolatedX, &x));
3559371c9d4SSatish 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])));
3569566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(interpolatedX, &x));
3579566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&interpolatedX));
358c4762a1bSJed Brown     user->next_output += PetscRealConstant(0.1);
359c4762a1bSJed Brown   }
3603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
361c4762a1bSJed Brown }
362c4762a1bSJed Brown 
main(int argc,char ** argv)363d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
364d71ae5a4SJacob Faibussowitsch {
365c4762a1bSJed Brown   Vec                P;
366c4762a1bSJed Brown   PetscBool          monitor = PETSC_FALSE;
367c4762a1bSJed Brown   PetscScalar       *x_ptr;
368c4762a1bSJed Brown   const PetscScalar *y_ptr;
369c4762a1bSJed Brown   PetscMPIInt        size;
370c4762a1bSJed Brown   struct _n_User     user;
371c4762a1bSJed Brown   Tao                tao;
372c4762a1bSJed Brown   KSP                ksp;
373c4762a1bSJed Brown   PC                 pc;
374c4762a1bSJed Brown 
375c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
376c4762a1bSJed Brown      Initialize program
377c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
378327415f7SBarry Smith   PetscFunctionBeginUser;
3799566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
3813c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
382c4762a1bSJed Brown 
383c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
3849566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
3859566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOBQNLS));
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
388c4762a1bSJed Brown     Set runtime options
389c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
390c4762a1bSJed Brown   user.next_output  = 0.0;
391c4762a1bSJed Brown   user.mu           = PetscRealConstant(1.0e3);
392c4762a1bSJed Brown   user.ftime        = PetscRealConstant(0.5);
393c4762a1bSJed Brown   user.implicitform = PETSC_TRUE;
394c4762a1bSJed Brown 
3959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
3969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
3979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &user.implicitform, NULL));
398c4762a1bSJed Brown 
399c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
4009566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.A));
4019566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
4029566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.A));
4039566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.A));
4049566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.U, NULL));
4059566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Lambda[0], NULL));
4069566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Lambda2[0], NULL));
4079566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Ihp1[0], NULL));
4089566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A, &user.Ihp2[0], NULL));
409c4762a1bSJed Brown 
4109566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Jacp));
4119566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
4129566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.Jacp));
4139566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.Jacp));
4149566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &user.Dir, NULL));
4159566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &user.Mup[0], NULL));
4169566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &user.Mup2[0], NULL));
4179566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &user.Ihp3[0], NULL));
4189566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &user.Ihp4[0], NULL));
419c4762a1bSJed Brown 
420c4762a1bSJed Brown   /* Create timestepping solver context */
4219566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &user.ts));
4229566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(user.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
423c4762a1bSJed Brown   if (user.implicitform) {
4249566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(user.ts, NULL, IFunction, &user));
4259566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(user.ts, user.A, user.A, IJacobian, &user));
4269566063dSJacob Faibussowitsch     PetscCall(TSSetType(user.ts, TSCN));
427c4762a1bSJed Brown   } else {
4289566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(user.ts, NULL, RHSFunction, &user));
4299566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(user.ts, user.A, user.A, RHSJacobian, &user));
4309566063dSJacob Faibussowitsch     PetscCall(TSSetType(user.ts, TSRK));
431c4762a1bSJed Brown   }
4329566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(user.ts, user.Jacp, RHSJacobianP, &user));
4339566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(user.ts, user.ftime));
4349566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(user.ts, TS_EXACTFINALTIME_MATCHSTEP));
435c4762a1bSJed Brown 
43648a46eb9SPierre Jolivet   if (monitor) PetscCall(TSMonitorSet(user.ts, Monitor, &user, NULL));
437c4762a1bSJed Brown 
438c4762a1bSJed Brown   /* Set ODE initial conditions */
4399566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.U, &x_ptr));
440c4762a1bSJed Brown   x_ptr[0] = 2.0;
441c4762a1bSJed Brown   x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu);
4429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.U, &x_ptr));
4439566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(user.ts, PetscRealConstant(0.001)));
444c4762a1bSJed Brown 
445c4762a1bSJed Brown   /* Set runtime options */
4469566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(user.ts));
447c4762a1bSJed Brown 
4489566063dSJacob Faibussowitsch   PetscCall(TSSolve(user.ts, user.U));
4499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user.U, &y_ptr));
450c4762a1bSJed Brown   user.ob[0] = y_ptr[0];
451c4762a1bSJed Brown   user.ob[1] = y_ptr[1];
4529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user.U, &y_ptr));
453c4762a1bSJed Brown 
454c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used.
455c4762a1bSJed Brown      Skip checkpointing for the first TSSolve since no adjoint run follows it.
456c4762a1bSJed Brown    */
4579566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(user.ts));
458c4762a1bSJed Brown 
459c4762a1bSJed Brown   /* Optimization starts */
4609566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.H));
4619566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.H, PETSC_DECIDE, PETSC_DECIDE, 1, 1));
4629566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.H)); /* Hessian should be symmetric. Do we need to do MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE) ? */
463c4762a1bSJed Brown 
464c4762a1bSJed Brown   /* Set initial solution guess */
4659566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &P, NULL));
4669566063dSJacob Faibussowitsch   PetscCall(VecGetArray(P, &x_ptr));
467c4762a1bSJed Brown   x_ptr[0] = PetscRealConstant(1.2);
4689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(P, &x_ptr));
4699566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, P));
470c4762a1bSJed Brown 
471c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
4729566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
4739566063dSJacob Faibussowitsch   PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user));
474c4762a1bSJed Brown 
475c4762a1bSJed Brown   /* Check for any TAO command line options */
4769566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao, &ksp));
477c4762a1bSJed Brown   if (ksp) {
4789566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
4799566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, PCNONE));
480c4762a1bSJed Brown   }
4819566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
482c4762a1bSJed Brown 
4839566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
484c4762a1bSJed Brown 
4859566063dSJacob Faibussowitsch   PetscCall(VecView(P, PETSC_VIEWER_STDOUT_WORLD));
486c4762a1bSJed Brown   /* Free TAO data structures */
4879566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
488c4762a1bSJed Brown 
489c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
490c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
491c4762a1bSJed Brown      are no longer needed.
492c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.H));
4949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.A));
4959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.U));
4969566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Jacp));
4979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Lambda[0]));
4989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Mup[0]));
4999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Lambda2[0]));
5009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Mup2[0]));
5019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp1[0]));
5029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp2[0]));
5039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp3[0]));
5049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp4[0]));
5059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Dir));
5069566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&user.ts));
5079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&P));
5089566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
509b122ec5aSJacob Faibussowitsch   return 0;
510c4762a1bSJed Brown }
511c4762a1bSJed Brown 
512c4762a1bSJed Brown /* ------------------------------------------------------------------ */
513c4762a1bSJed Brown /*
514c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
515c4762a1bSJed Brown 
516c4762a1bSJed Brown    Input Parameters:
517c4762a1bSJed Brown    tao - the Tao context
518c4762a1bSJed Brown    X   - the input vector
519a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
520c4762a1bSJed Brown 
521c4762a1bSJed Brown    Output Parameters:
522c4762a1bSJed Brown    f   - the newly evaluated function
523c4762a1bSJed Brown    G   - the newly evaluated gradient
524c4762a1bSJed Brown */
FormFunctionGradient(Tao tao,Vec P,PetscReal * f,Vec G,PetscCtx ctx)525*2a8381b2SBarry Smith PetscErrorCode FormFunctionGradient(Tao tao, Vec P, PetscReal *f, Vec G, PetscCtx ctx)
526d71ae5a4SJacob Faibussowitsch {
527c4762a1bSJed Brown   User               user_ptr = (User)ctx;
528c4762a1bSJed Brown   TS                 ts       = user_ptr->ts;
529c4762a1bSJed Brown   PetscScalar       *x_ptr, *g;
530c4762a1bSJed Brown   const PetscScalar *y_ptr;
531c4762a1bSJed Brown 
532c4762a1bSJed Brown   PetscFunctionBeginUser;
5339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P, &y_ptr));
534c4762a1bSJed Brown   user_ptr->mu = y_ptr[0];
5359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P, &y_ptr));
536c4762a1bSJed Brown 
5379566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
5389566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
5399566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, PetscRealConstant(0.001))); /* can be overwritten by command line options */
5409566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
5419566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->U, &x_ptr));
542c4762a1bSJed Brown   x_ptr[0] = 2.0;
543c4762a1bSJed 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);
5449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->U, &x_ptr));
545c4762a1bSJed Brown 
5469566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, user_ptr->U));
547c4762a1bSJed Brown 
5489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user_ptr->U, &y_ptr));
549c4762a1bSJed 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]);
550c4762a1bSJed Brown 
551c4762a1bSJed Brown   /*   Reset initial conditions for the adjoint integration */
5529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Lambda[0], &x_ptr));
553c4762a1bSJed Brown   x_ptr[0] = 2. * (y_ptr[0] - user_ptr->ob[0]);
554c4762a1bSJed Brown   x_ptr[1] = 2. * (y_ptr[1] - user_ptr->ob[1]);
5559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user_ptr->U, &y_ptr));
5569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Lambda[0], &x_ptr));
557c4762a1bSJed Brown 
5589566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Mup[0], &x_ptr));
559c4762a1bSJed Brown   x_ptr[0] = 0.0;
5609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Mup[0], &x_ptr));
5619566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, user_ptr->Lambda, user_ptr->Mup));
562c4762a1bSJed Brown 
5639566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
564c4762a1bSJed Brown 
5659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Mup[0], &x_ptr));
5669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user_ptr->Lambda[0], &y_ptr));
5679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G, &g));
568c4762a1bSJed 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];
5699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Mup[0], &x_ptr));
5709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user_ptr->Lambda[0], &y_ptr));
5719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G, &g));
5723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
573c4762a1bSJed Brown }
574c4762a1bSJed Brown 
FormHessian(Tao tao,Vec P,Mat H,Mat Hpre,PetscCtx ctx)575*2a8381b2SBarry Smith PetscErrorCode FormHessian(Tao tao, Vec P, Mat H, Mat Hpre, PetscCtx ctx)
576d71ae5a4SJacob Faibussowitsch {
577c4762a1bSJed Brown   User           user_ptr = (User)ctx;
578c4762a1bSJed Brown   PetscScalar    harr[1];
579c4762a1bSJed Brown   const PetscInt rows[1] = {0};
580c4762a1bSJed Brown   PetscInt       col     = 0;
581c4762a1bSJed Brown 
582c4762a1bSJed Brown   PetscFunctionBeginUser;
5839566063dSJacob Faibussowitsch   PetscCall(Adjoint2(P, harr, user_ptr));
5849566063dSJacob Faibussowitsch   PetscCall(MatSetValues(H, 1, rows, 1, &col, harr, INSERT_VALUES));
585c4762a1bSJed Brown 
5869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
5879566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
588c4762a1bSJed Brown   if (H != Hpre) {
5899566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Hpre, MAT_FINAL_ASSEMBLY));
5909566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Hpre, MAT_FINAL_ASSEMBLY));
591c4762a1bSJed Brown   }
5923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
593c4762a1bSJed Brown }
594c4762a1bSJed Brown 
Adjoint2(Vec P,PetscScalar arr[],User ctx)595d71ae5a4SJacob Faibussowitsch PetscErrorCode Adjoint2(Vec P, PetscScalar arr[], User ctx)
596d71ae5a4SJacob Faibussowitsch {
597c4762a1bSJed Brown   TS                 ts = ctx->ts;
598c4762a1bSJed Brown   const PetscScalar *z_ptr;
599c4762a1bSJed Brown   PetscScalar       *x_ptr, *y_ptr, dzdp, dzdp2;
600c4762a1bSJed Brown   Mat                tlmsen;
601c4762a1bSJed Brown 
602c4762a1bSJed Brown   PetscFunctionBeginUser;
603c4762a1bSJed Brown   /* Reset TSAdjoint so that AdjointSetUp will be called again */
6049566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
605c4762a1bSJed Brown 
606c4762a1bSJed Brown   /* The directional vector should be 1 since it is one-dimensional */
6079566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Dir, &x_ptr));
608c4762a1bSJed Brown   x_ptr[0] = 1.;
6099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Dir, &x_ptr));
610c4762a1bSJed Brown 
6119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P, &z_ptr));
612c4762a1bSJed Brown   ctx->mu = z_ptr[0];
6139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P, &z_ptr));
614c4762a1bSJed Brown 
615c4762a1bSJed Brown   dzdp  = -10.0 / (81.0 * ctx->mu * ctx->mu) + 2.0 * 292.0 / (2187.0 * ctx->mu * ctx->mu * ctx->mu);
616c4762a1bSJed 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);
617c4762a1bSJed Brown 
6189566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
6199566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
6209566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, PetscRealConstant(0.001))); /* can be overwritten by command line options */
6219566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
6229566063dSJacob Faibussowitsch   PetscCall(TSSetCostHessianProducts(ts, 1, ctx->Lambda2, ctx->Mup2, ctx->Dir));
623c4762a1bSJed Brown 
6249566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(ctx->Jacp));
6259566063dSJacob Faibussowitsch   PetscCall(MatSetValue(ctx->Jacp, 1, 0, dzdp, INSERT_VALUES));
6269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(ctx->Jacp, MAT_FINAL_ASSEMBLY));
6279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(ctx->Jacp, MAT_FINAL_ASSEMBLY));
628c4762a1bSJed Brown 
6299566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetForward(ts, ctx->Jacp));
6309566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->U, &y_ptr));
631c4762a1bSJed Brown   y_ptr[0] = 2.0;
632c4762a1bSJed Brown   y_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * ctx->mu) - 292.0 / (2187.0 * ctx->mu * ctx->mu);
6339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->U, &y_ptr));
6349566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, ctx->U));
635c4762a1bSJed Brown 
636c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
6379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->U, &z_ptr));
6389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda[0], &y_ptr));
639c4762a1bSJed Brown   y_ptr[0] = 2. * (z_ptr[0] - ctx->ob[0]);
640c4762a1bSJed Brown   y_ptr[1] = 2. * (z_ptr[1] - ctx->ob[1]);
6419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda[0], &y_ptr));
6429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->U, &z_ptr));
6439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Mup[0], &y_ptr));
644c4762a1bSJed Brown   y_ptr[0] = 0.0;
6459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Mup[0], &y_ptr));
6469566063dSJacob Faibussowitsch   PetscCall(TSForwardGetSensitivities(ts, NULL, &tlmsen));
6479566063dSJacob Faibussowitsch   PetscCall(MatDenseGetColumn(tlmsen, 0, &x_ptr));
6489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda2[0], &y_ptr));
649c4762a1bSJed Brown   y_ptr[0] = 2. * x_ptr[0];
650c4762a1bSJed Brown   y_ptr[1] = 2. * x_ptr[1];
6519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0], &y_ptr));
6529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Mup2[0], &y_ptr));
653c4762a1bSJed Brown   y_ptr[0] = 0.0;
6549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Mup2[0], &y_ptr));
6559566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreColumn(tlmsen, &x_ptr));
6569566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, ctx->Lambda, ctx->Mup));
657c4762a1bSJed Brown   if (ctx->implicitform) {
6589566063dSJacob Faibussowitsch     PetscCall(TSSetIHessianProduct(ts, ctx->Ihp1, IHessianProductUU, ctx->Ihp2, IHessianProductUP, ctx->Ihp3, IHessianProductPU, ctx->Ihp4, IHessianProductPP, ctx));
659c4762a1bSJed Brown   } else {
6609566063dSJacob Faibussowitsch     PetscCall(TSSetRHSHessianProduct(ts, ctx->Ihp1, RHSHessianProductUU, ctx->Ihp2, RHSHessianProductUP, ctx->Ihp3, RHSHessianProductPU, ctx->Ihp4, RHSHessianProductPP, ctx));
661c4762a1bSJed Brown   }
6629566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
663c4762a1bSJed Brown 
6649566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda[0], &x_ptr));
6659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda2[0], &y_ptr));
6669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->Mup2[0], &z_ptr));
667c4762a1bSJed Brown 
668c4762a1bSJed Brown   arr[0] = x_ptr[1] * dzdp2 + y_ptr[1] * dzdp2 + z_ptr[0];
669c4762a1bSJed Brown 
6709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0], &x_ptr));
6719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0], &y_ptr));
6729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->Mup2[0], &z_ptr));
673c4762a1bSJed Brown 
674c4762a1bSJed Brown   /* Disable second-order adjoint mode */
6759566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
6769566063dSJacob Faibussowitsch   PetscCall(TSAdjointResetForward(ts));
6773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
678c4762a1bSJed Brown }
679c4762a1bSJed Brown 
680c4762a1bSJed Brown /*TEST
681c4762a1bSJed Brown     build:
682c4762a1bSJed Brown       requires: !complex !single
683c4762a1bSJed Brown     test:
684188af4bfSBarry Smith       args: -implicitform 0 -ts_type rk -ts_adapt_type none -mu 10 -ts_time_step 0.1 -viewer_binary_skip_info -tao_monitor -tao_view
685c4762a1bSJed Brown       output_file: output/ex20opt_p_1.out
686c4762a1bSJed Brown 
687c4762a1bSJed Brown     test:
688c4762a1bSJed Brown       suffix: 2
689188af4bfSBarry Smith       args: -implicitform 0 -ts_type rk -ts_adapt_type none -mu 10 -ts_time_step 0.01 -viewer_binary_skip_info -tao_monitor -tao_type bntr -tao_bnk_pc_type none
690c4762a1bSJed Brown       output_file: output/ex20opt_p_2.out
691c4762a1bSJed Brown 
692c4762a1bSJed Brown     test:
693c4762a1bSJed Brown       suffix: 3
694188af4bfSBarry Smith       args: -ts_type cn -ts_adapt_type none -mu 100 -ts_time_step 0.01 -viewer_binary_skip_info -tao_monitor -tao_view
695c4762a1bSJed Brown       output_file: output/ex20opt_p_3.out
696c4762a1bSJed Brown 
697c4762a1bSJed Brown     test:
698c4762a1bSJed Brown       suffix: 4
699188af4bfSBarry Smith       args: -ts_type cn -ts_adapt_type none -mu 100 -ts_time_step 0.01 -viewer_binary_skip_info -tao_monitor -tao_type bntr -tao_bnk_pc_type none
700c4762a1bSJed Brown       output_file: output/ex20opt_p_4.out
701c4762a1bSJed Brown 
702c4762a1bSJed Brown TEST*/
703