xref: /petsc/src/ts/tutorials/ex20opt_p.c (revision 63a3b9bc7a1f24f247904ccba9383635fe6abade)
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 
44c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
45c4762a1bSJed Brown {
46c4762a1bSJed Brown   User              user = (User)ctx;
47c4762a1bSJed Brown   PetscScalar       *f;
48c4762a1bSJed Brown   const PetscScalar *u;
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   PetscFunctionBeginUser;
519566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
53c4762a1bSJed Brown   f[0] = u[1];
54c4762a1bSJed Brown   f[1] = user->mu*((1.-u[0]*u[0])*u[1]-u[0]);
559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
57c4762a1bSJed Brown   PetscFunctionReturn(0);
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
60c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
61c4762a1bSJed Brown {
62c4762a1bSJed Brown   User              user = (User)ctx;
63c4762a1bSJed Brown   PetscReal         mu   = user->mu;
64c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
65c4762a1bSJed Brown   PetscScalar       J[2][2];
66c4762a1bSJed Brown   const PetscScalar *u;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   PetscFunctionBeginUser;
699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   J[0][0] = 0;
72c4762a1bSJed Brown   J[1][0] = -mu*(2.0*u[1]*u[0]+1.);
73c4762a1bSJed Brown   J[0][1] = 1.0;
74c4762a1bSJed Brown   J[1][1] = mu*(1.0-u[0]*u[0]);
759566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
78c4762a1bSJed Brown   if (B && A != B) {
799566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
809566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
81c4762a1bSJed Brown   }
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
84c4762a1bSJed Brown   PetscFunctionReturn(0);
85c4762a1bSJed Brown }
86c4762a1bSJed Brown 
87c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
88c4762a1bSJed Brown {
89c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
90c4762a1bSJed Brown   PetscScalar       *vhv;
91c4762a1bSJed Brown   PetscScalar       dJdU[2][2][2]={{{0}}};
92c4762a1bSJed Brown   PetscInt          i,j,k;
93c4762a1bSJed Brown   User              user = (User)ctx;
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   PetscFunctionBeginUser;
969566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
979566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0],&vl));
989566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr,&vr));
999566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0],&vhv));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   dJdU[1][0][0] = -2.*user->mu*u[1];
102c4762a1bSJed Brown   dJdU[1][1][0] = -2.*user->mu*u[0];
103c4762a1bSJed Brown   dJdU[1][0][1] = -2.*user->mu*u[0];
104c4762a1bSJed Brown   for (j=0; j<2; j++) {
105c4762a1bSJed Brown     vhv[j] = 0;
106c4762a1bSJed Brown     for (k=0; k<2; k++)
107c4762a1bSJed Brown       for (i=0; i<2; i++)
108c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
109c4762a1bSJed Brown   }
110c4762a1bSJed Brown 
1119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0],&vl));
1139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr,&vr));
1149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0],&vhv));
115c4762a1bSJed Brown   PetscFunctionReturn(0);
116c4762a1bSJed Brown }
117c4762a1bSJed Brown 
118c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
119c4762a1bSJed Brown {
120c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
121c4762a1bSJed Brown   PetscScalar       *vhv;
122c4762a1bSJed Brown   PetscScalar       dJdP[2][2][1]={{{0}}};
123c4762a1bSJed Brown   PetscInt          i,j,k;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   PetscFunctionBeginUser;
1269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0],&vl));
1289566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr,&vr));
1299566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0],&vhv));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   dJdP[1][0][0] = -(1.+2.*u[0]*u[1]);
132c4762a1bSJed Brown   dJdP[1][1][0] = 1.-u[0]*u[0];
133c4762a1bSJed Brown   for (j=0; j<2; j++) {
134c4762a1bSJed Brown     vhv[j] = 0;
135c4762a1bSJed Brown     for (k=0; k<1; k++)
136c4762a1bSJed Brown       for (i=0; i<2; i++)
137c4762a1bSJed Brown         vhv[j] += vl[i]*dJdP[i][j][k]*vr[k];
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
1409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0],&vl));
1429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr,&vr));
1439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0],&vhv));
144c4762a1bSJed Brown   PetscFunctionReturn(0);
145c4762a1bSJed Brown }
146c4762a1bSJed Brown 
147c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
148c4762a1bSJed Brown {
149c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
150c4762a1bSJed Brown   PetscScalar       *vhv;
151c4762a1bSJed Brown   PetscScalar       dJdU[2][1][2]={{{0}}};
152c4762a1bSJed Brown   PetscInt          i,j,k;
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   PetscFunctionBeginUser;
1559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0],&vl));
1579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr,&vr));
1589566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0],&vhv));
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   dJdU[1][0][0] = -1.-2.*u[1]*u[0];
161c4762a1bSJed Brown   dJdU[1][0][1] = 1.-u[0]*u[0];
162c4762a1bSJed Brown   for (j=0; j<1; j++) {
163c4762a1bSJed Brown     vhv[j] = 0;
164c4762a1bSJed Brown     for (k=0; k<2; k++)
165c4762a1bSJed Brown       for (i=0; i<2; i++)
166c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
167c4762a1bSJed Brown   }
168c4762a1bSJed Brown 
1699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0],&vl));
1719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr,&vr));
1729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0],&vhv));
173c4762a1bSJed Brown   PetscFunctionReturn(0);
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
177c4762a1bSJed Brown {
178c4762a1bSJed Brown   PetscFunctionBeginUser;
179c4762a1bSJed Brown   PetscFunctionReturn(0);
180c4762a1bSJed Brown }
181c4762a1bSJed Brown 
182c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE  -------------------- */
183c4762a1bSJed Brown 
184c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
185c4762a1bSJed Brown {
186c4762a1bSJed Brown   User              user = (User)ctx;
187c4762a1bSJed Brown   PetscScalar       *f;
188c4762a1bSJed Brown   const PetscScalar *u,*udot;
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   PetscFunctionBeginUser;
1919566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1929566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
1939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   f[0] = udot[0] - u[1];
196c4762a1bSJed Brown   f[1] = udot[1] - user->mu*((1.0-u[0]*u[0])*u[1] - u[0]) ;
197c4762a1bSJed Brown 
1989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot,&udot));
2009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
201c4762a1bSJed Brown   PetscFunctionReturn(0);
202c4762a1bSJed Brown }
203c4762a1bSJed Brown 
204c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
205c4762a1bSJed Brown {
206c4762a1bSJed Brown   User              user     = (User)ctx;
207c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
208c4762a1bSJed Brown   PetscScalar       J[2][2];
209c4762a1bSJed Brown   const PetscScalar *u;
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   PetscFunctionBeginUser;
2129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   J[0][0] = a;     J[0][1] =  -1.0;
215c4762a1bSJed Brown   J[1][0] = user->mu*(1.0 + 2.0*u[0]*u[1]);   J[1][1] = a - user->mu*(1.0-u[0]*u[0]);
2169566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
2179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
2189566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
2199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
220c4762a1bSJed Brown   if (A != B) {
2219566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
2229566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
223c4762a1bSJed Brown   }
224c4762a1bSJed Brown 
2259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
226c4762a1bSJed Brown   PetscFunctionReturn(0);
227c4762a1bSJed Brown }
228c4762a1bSJed Brown 
229c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec U,Mat A,void *ctx)
230c4762a1bSJed Brown {
231c4762a1bSJed Brown   PetscInt          row[] = {0,1},col[]={0};
232c4762a1bSJed Brown   PetscScalar       J[2][1];
233c4762a1bSJed Brown   const PetscScalar *u;
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   PetscFunctionBeginUser;
2369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
237c4762a1bSJed Brown 
238c4762a1bSJed Brown   J[0][0] = 0;
239c4762a1bSJed Brown   J[1][0] = (1.-u[0]*u[0])*u[1]-u[0];
2409566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES));
2419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
2429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
2439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
244c4762a1bSJed Brown 
2459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
246c4762a1bSJed Brown   PetscFunctionReturn(0);
247c4762a1bSJed Brown }
248c4762a1bSJed Brown 
249c4762a1bSJed Brown static PetscErrorCode IHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
250c4762a1bSJed Brown {
251c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
252c4762a1bSJed Brown   PetscScalar       *vhv;
253c4762a1bSJed Brown   PetscScalar       dJdU[2][2][2]={{{0}}};
254c4762a1bSJed Brown   PetscInt          i,j,k;
255c4762a1bSJed Brown   User              user = (User)ctx;
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   PetscFunctionBeginUser;
2589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
2599566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0],&vl));
2609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr,&vr));
2619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0],&vhv));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   dJdU[1][0][0] = 2.*user->mu*u[1];
264c4762a1bSJed Brown   dJdU[1][1][0] = 2.*user->mu*u[0];
265c4762a1bSJed Brown   dJdU[1][0][1] = 2.*user->mu*u[0];
266c4762a1bSJed Brown   for (j=0; j<2; j++) {
267c4762a1bSJed Brown     vhv[j] = 0;
268c4762a1bSJed Brown     for (k=0; k<2; k++)
269c4762a1bSJed Brown       for (i=0; i<2; i++)
270c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
271c4762a1bSJed Brown   }
272c4762a1bSJed Brown 
2739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
2749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0],&vl));
2759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr,&vr));
2769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0],&vhv));
277c4762a1bSJed Brown   PetscFunctionReturn(0);
278c4762a1bSJed Brown }
279c4762a1bSJed Brown 
280c4762a1bSJed Brown static PetscErrorCode IHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
281c4762a1bSJed Brown {
282c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
283c4762a1bSJed Brown   PetscScalar       *vhv;
284c4762a1bSJed Brown   PetscScalar       dJdP[2][2][1]={{{0}}};
285c4762a1bSJed Brown   PetscInt          i,j,k;
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   PetscFunctionBeginUser;
2889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
2899566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0],&vl));
2909566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr,&vr));
2919566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0],&vhv));
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   dJdP[1][0][0] = 1.+2.*u[0]*u[1];
294c4762a1bSJed Brown   dJdP[1][1][0] = u[0]*u[0]-1.;
295c4762a1bSJed Brown   for (j=0; j<2; j++) {
296c4762a1bSJed Brown     vhv[j] = 0;
297c4762a1bSJed Brown     for (k=0; k<1; k++)
298c4762a1bSJed Brown       for (i=0; i<2; i++)
299c4762a1bSJed Brown         vhv[j] += vl[i]*dJdP[i][j][k]*vr[k];
300c4762a1bSJed Brown   }
301c4762a1bSJed Brown 
3029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
3039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0],&vl));
3049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr,&vr));
3059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0],&vhv));
306c4762a1bSJed Brown   PetscFunctionReturn(0);
307c4762a1bSJed Brown }
308c4762a1bSJed Brown 
309c4762a1bSJed Brown static PetscErrorCode IHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
310c4762a1bSJed Brown {
311c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
312c4762a1bSJed Brown   PetscScalar       *vhv;
313c4762a1bSJed Brown   PetscScalar       dJdU[2][1][2]={{{0}}};
314c4762a1bSJed Brown   PetscInt          i,j,k;
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   PetscFunctionBeginUser;
3179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
3189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0],&vl));
3199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr,&vr));
3209566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0],&vhv));
321c4762a1bSJed Brown 
322c4762a1bSJed Brown   dJdU[1][0][0] = 1.+2.*u[1]*u[0];
323c4762a1bSJed Brown   dJdU[1][0][1] = u[0]*u[0]-1.;
324c4762a1bSJed Brown   for (j=0; j<1; j++) {
325c4762a1bSJed Brown     vhv[j] = 0;
326c4762a1bSJed Brown     for (k=0; k<2; k++)
327c4762a1bSJed Brown       for (i=0; i<2; i++)
328c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
329c4762a1bSJed Brown   }
330c4762a1bSJed Brown 
3319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
3329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0],&vl));
3339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr,&vr));
3349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0],&vhv));
335c4762a1bSJed Brown   PetscFunctionReturn(0);
336c4762a1bSJed Brown }
337c4762a1bSJed Brown 
338c4762a1bSJed Brown static PetscErrorCode IHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
339c4762a1bSJed Brown {
340c4762a1bSJed Brown   PetscFunctionBeginUser;
341c4762a1bSJed Brown   PetscFunctionReturn(0);
342c4762a1bSJed Brown }
343c4762a1bSJed Brown 
344c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
345c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
346c4762a1bSJed Brown {
347c4762a1bSJed Brown   const PetscScalar *x;
348c4762a1bSJed Brown   PetscReal         tfinal, dt;
349c4762a1bSJed Brown   User              user = (User)ctx;
350c4762a1bSJed Brown   Vec               interpolatedX;
351c4762a1bSJed Brown 
352c4762a1bSJed Brown   PetscFunctionBeginUser;
3539566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts,&dt));
3549566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts,&tfinal));
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
3579566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X,&interpolatedX));
3589566063dSJacob Faibussowitsch     PetscCall(TSInterpolate(ts,user->next_output,interpolatedX));
3599566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(interpolatedX,&x));
360*63a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[%g] %" PetscInt_FMT " TS %g (dt = %g) X %g %g\n",
361*63a3b9bcSJacob Faibussowitsch                           (double)user->next_output,step,(double)t,(double)dt,
362*63a3b9bcSJacob Faibussowitsch                           (double)PetscRealPart(x[0]),(double)PetscRealPart(x[1])));
3639566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(interpolatedX,&x));
3649566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&interpolatedX));
365c4762a1bSJed Brown     user->next_output += PetscRealConstant(0.1);
366c4762a1bSJed Brown   }
367c4762a1bSJed Brown   PetscFunctionReturn(0);
368c4762a1bSJed Brown }
369c4762a1bSJed Brown 
370c4762a1bSJed Brown int main(int argc,char **argv)
371c4762a1bSJed Brown {
372c4762a1bSJed Brown   Vec                P;
373c4762a1bSJed Brown   PetscBool          monitor = PETSC_FALSE;
374c4762a1bSJed Brown   PetscScalar        *x_ptr;
375c4762a1bSJed Brown   const PetscScalar  *y_ptr;
376c4762a1bSJed Brown   PetscMPIInt        size;
377c4762a1bSJed Brown   struct _n_User     user;
378c4762a1bSJed Brown   Tao                tao;
379c4762a1bSJed Brown   KSP                ksp;
380c4762a1bSJed Brown   PC                 pc;
381c4762a1bSJed Brown 
382c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
383c4762a1bSJed Brown      Initialize program
384c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3859566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
3869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
3873c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
388c4762a1bSJed Brown 
389c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
3909566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
3919566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOBQNLS));
392c4762a1bSJed Brown 
393c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
394c4762a1bSJed Brown     Set runtime options
395c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
396c4762a1bSJed Brown   user.next_output  = 0.0;
397c4762a1bSJed Brown   user.mu           = PetscRealConstant(1.0e3);
398c4762a1bSJed Brown   user.ftime        = PetscRealConstant(0.5);
399c4762a1bSJed Brown   user.implicitform = PETSC_TRUE;
400c4762a1bSJed Brown 
4019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
4029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
4039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-implicitform",&user.implicitform,NULL));
404c4762a1bSJed Brown 
405c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
4069566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user.A));
4079566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2));
4089566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.A));
4099566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.A));
4109566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A,&user.U,NULL));
4119566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A,&user.Lambda[0],NULL));
4129566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A,&user.Lambda2[0],NULL));
4139566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A,&user.Ihp1[0],NULL));
4149566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A,&user.Ihp2[0],NULL));
415c4762a1bSJed Brown 
4169566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Jacp));
4179566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1));
4189566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.Jacp));
4199566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.Jacp));
4209566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp,&user.Dir,NULL));
4219566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp,&user.Mup[0],NULL));
4229566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp,&user.Mup2[0],NULL));
4239566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp,&user.Ihp3[0],NULL));
4249566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp,&user.Ihp4[0],NULL));
425c4762a1bSJed Brown 
426c4762a1bSJed Brown   /* Create timestepping solver context */
4279566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&user.ts));
4289566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(user.ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
429c4762a1bSJed Brown   if (user.implicitform) {
4309566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(user.ts,NULL,IFunction,&user));
4319566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(user.ts,user.A,user.A,IJacobian,&user));
4329566063dSJacob Faibussowitsch     PetscCall(TSSetType(user.ts,TSCN));
433c4762a1bSJed Brown   } else {
4349566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(user.ts,NULL,RHSFunction,&user));
4359566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(user.ts,user.A,user.A,RHSJacobian,&user));
4369566063dSJacob Faibussowitsch     PetscCall(TSSetType(user.ts,TSRK));
437c4762a1bSJed Brown   }
4389566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(user.ts,user.Jacp,RHSJacobianP,&user));
4399566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(user.ts,user.ftime));
4409566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(user.ts,TS_EXACTFINALTIME_MATCHSTEP));
441c4762a1bSJed Brown 
442c4762a1bSJed Brown   if (monitor) {
4439566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(user.ts,Monitor,&user,NULL));
444c4762a1bSJed Brown   }
445c4762a1bSJed Brown 
446c4762a1bSJed Brown   /* Set ODE initial conditions */
4479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.U,&x_ptr));
448c4762a1bSJed Brown   x_ptr[0] = 2.0;
449c4762a1bSJed Brown   x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu);
4509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.U,&x_ptr));
4519566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(user.ts,PetscRealConstant(0.001)));
452c4762a1bSJed Brown 
453c4762a1bSJed Brown   /* Set runtime options */
4549566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(user.ts));
455c4762a1bSJed Brown 
4569566063dSJacob Faibussowitsch   PetscCall(TSSolve(user.ts,user.U));
4579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user.U,&y_ptr));
458c4762a1bSJed Brown   user.ob[0] = y_ptr[0];
459c4762a1bSJed Brown   user.ob[1] = y_ptr[1];
4609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user.U,&y_ptr));
461c4762a1bSJed Brown 
462c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used.
463c4762a1bSJed Brown      Skip checkpointing for the first TSSolve since no adjoint run follows it.
464c4762a1bSJed Brown    */
4659566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(user.ts));
466c4762a1bSJed Brown 
467c4762a1bSJed Brown   /* Optimization starts */
4689566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user.H));
4699566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.H,PETSC_DECIDE,PETSC_DECIDE,1,1));
4709566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.H)); /* Hessian should be symmetric. Do we need to do MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE) ? */
471c4762a1bSJed Brown 
472c4762a1bSJed Brown   /* Set initial solution guess */
4739566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp,&P,NULL));
4749566063dSJacob Faibussowitsch   PetscCall(VecGetArray(P,&x_ptr));
475c4762a1bSJed Brown   x_ptr[0] = PetscRealConstant(1.2);
4769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(P,&x_ptr));
4779566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,P));
478c4762a1bSJed Brown 
479c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
4809566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user));
4819566063dSJacob Faibussowitsch   PetscCall(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user));
482c4762a1bSJed Brown 
483c4762a1bSJed Brown   /* Check for any TAO command line options */
4849566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao,&ksp));
485c4762a1bSJed Brown   if (ksp) {
4869566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp,&pc));
4879566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc,PCNONE));
488c4762a1bSJed Brown   }
4899566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
490c4762a1bSJed Brown 
4919566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
492c4762a1bSJed Brown 
4939566063dSJacob Faibussowitsch   PetscCall(VecView(P,PETSC_VIEWER_STDOUT_WORLD));
494c4762a1bSJed Brown   /* Free TAO data structures */
4959566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
496c4762a1bSJed Brown 
497c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
498c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
499c4762a1bSJed Brown      are no longer needed.
500c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
5019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.H));
5029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.A));
5039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.U));
5049566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Jacp));
5059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Lambda[0]));
5069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Mup[0]));
5079566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Lambda2[0]));
5089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Mup2[0]));
5099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp1[0]));
5109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp2[0]));
5119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp3[0]));
5129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Ihp4[0]));
5139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.Dir));
5149566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&user.ts));
5159566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&P));
5169566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
517b122ec5aSJacob Faibussowitsch   return 0;
518c4762a1bSJed Brown }
519c4762a1bSJed Brown 
520c4762a1bSJed Brown /* ------------------------------------------------------------------ */
521c4762a1bSJed Brown /*
522c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
523c4762a1bSJed Brown 
524c4762a1bSJed Brown    Input Parameters:
525c4762a1bSJed Brown    tao - the Tao context
526c4762a1bSJed Brown    X   - the input vector
527a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
528c4762a1bSJed Brown 
529c4762a1bSJed Brown    Output Parameters:
530c4762a1bSJed Brown    f   - the newly evaluated function
531c4762a1bSJed Brown    G   - the newly evaluated gradient
532c4762a1bSJed Brown */
533c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
534c4762a1bSJed Brown {
535c4762a1bSJed Brown   User              user_ptr = (User)ctx;
536c4762a1bSJed Brown   TS                ts = user_ptr->ts;
537c4762a1bSJed Brown   PetscScalar       *x_ptr,*g;
538c4762a1bSJed Brown   const PetscScalar *y_ptr;
539c4762a1bSJed Brown 
540c4762a1bSJed Brown   PetscFunctionBeginUser;
5419566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P,&y_ptr));
542c4762a1bSJed Brown   user_ptr->mu = y_ptr[0];
5439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P,&y_ptr));
544c4762a1bSJed Brown 
5459566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts,0.0));
5469566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts,0));
5479566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,PetscRealConstant(0.001))); /* can be overwritten by command line options */
5489566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
5499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->U,&x_ptr));
550c4762a1bSJed Brown   x_ptr[0] = 2.0;
551c4762a1bSJed 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);
5529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->U,&x_ptr));
553c4762a1bSJed Brown 
5549566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,user_ptr->U));
555c4762a1bSJed Brown 
5569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user_ptr->U,&y_ptr));
557c4762a1bSJed 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]);
558c4762a1bSJed Brown 
559c4762a1bSJed Brown   /*   Reset initial conditions for the adjoint integration */
5609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Lambda[0],&x_ptr));
561c4762a1bSJed Brown   x_ptr[0] = 2.*(y_ptr[0]-user_ptr->ob[0]);
562c4762a1bSJed Brown   x_ptr[1] = 2.*(y_ptr[1]-user_ptr->ob[1]);
5639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user_ptr->U,&y_ptr));
5649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Lambda[0],&x_ptr));
565c4762a1bSJed Brown 
5669566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Mup[0],&x_ptr));
567c4762a1bSJed Brown   x_ptr[0] = 0.0;
5689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Mup[0],&x_ptr));
5699566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts,1,user_ptr->Lambda,user_ptr->Mup));
570c4762a1bSJed Brown 
5719566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
572c4762a1bSJed Brown 
5739566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user_ptr->Mup[0],&x_ptr));
5749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user_ptr->Lambda[0],&y_ptr));
5759566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G,&g));
576c4762a1bSJed 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];
5779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user_ptr->Mup[0],&x_ptr));
5789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user_ptr->Lambda[0],&y_ptr));
5799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G,&g));
580c4762a1bSJed Brown   PetscFunctionReturn(0);
581c4762a1bSJed Brown }
582c4762a1bSJed Brown 
583c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec P,Mat H,Mat Hpre,void *ctx)
584c4762a1bSJed Brown {
585c4762a1bSJed Brown   User           user_ptr = (User)ctx;
586c4762a1bSJed Brown   PetscScalar    harr[1];
587c4762a1bSJed Brown   const PetscInt rows[1] = {0};
588c4762a1bSJed Brown   PetscInt       col = 0;
589c4762a1bSJed Brown 
590c4762a1bSJed Brown   PetscFunctionBeginUser;
5919566063dSJacob Faibussowitsch   PetscCall(Adjoint2(P,harr,user_ptr));
5929566063dSJacob Faibussowitsch   PetscCall(MatSetValues(H,1,rows,1,&col,harr,INSERT_VALUES));
593c4762a1bSJed Brown 
5949566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
5959566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
596c4762a1bSJed Brown   if (H != Hpre) {
5979566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(Hpre,MAT_FINAL_ASSEMBLY));
5989566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(Hpre,MAT_FINAL_ASSEMBLY));
599c4762a1bSJed Brown   }
600c4762a1bSJed Brown   PetscFunctionReturn(0);
601c4762a1bSJed Brown }
602c4762a1bSJed Brown 
603c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec P,PetscScalar arr[],User ctx)
604c4762a1bSJed Brown {
605c4762a1bSJed Brown   TS                ts = ctx->ts;
606c4762a1bSJed Brown   const PetscScalar *z_ptr;
607c4762a1bSJed Brown   PetscScalar       *x_ptr,*y_ptr,dzdp,dzdp2;
608c4762a1bSJed Brown   Mat               tlmsen;
609c4762a1bSJed Brown 
610c4762a1bSJed Brown   PetscFunctionBeginUser;
611c4762a1bSJed Brown   /* Reset TSAdjoint so that AdjointSetUp will be called again */
6129566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
613c4762a1bSJed Brown 
614c4762a1bSJed Brown   /* The directional vector should be 1 since it is one-dimensional */
6159566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Dir,&x_ptr));
616c4762a1bSJed Brown   x_ptr[0] = 1.;
6179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Dir,&x_ptr));
618c4762a1bSJed Brown 
6199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P,&z_ptr));
620c4762a1bSJed Brown   ctx->mu = z_ptr[0];
6219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P,&z_ptr));
622c4762a1bSJed Brown 
623c4762a1bSJed Brown   dzdp  = -10.0/(81.0*ctx->mu*ctx->mu) + 2.0*292.0/(2187.0*ctx->mu*ctx->mu*ctx->mu);
624c4762a1bSJed 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);
625c4762a1bSJed Brown 
6269566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts,0.0));
6279566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts,0));
6289566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,PetscRealConstant(0.001))); /* can be overwritten by command line options */
6299566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
6309566063dSJacob Faibussowitsch   PetscCall(TSSetCostHessianProducts(ts,1,ctx->Lambda2,ctx->Mup2,ctx->Dir));
631c4762a1bSJed Brown 
6329566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(ctx->Jacp));
6339566063dSJacob Faibussowitsch   PetscCall(MatSetValue(ctx->Jacp,1,0,dzdp,INSERT_VALUES));
6349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(ctx->Jacp,MAT_FINAL_ASSEMBLY));
6359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(ctx->Jacp,MAT_FINAL_ASSEMBLY));
636c4762a1bSJed Brown 
6379566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetForward(ts,ctx->Jacp));
6389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->U,&y_ptr));
639c4762a1bSJed Brown   y_ptr[0] = 2.0;
640c4762a1bSJed Brown   y_ptr[1] = -2.0/3.0 + 10.0/(81.0*ctx->mu) - 292.0/(2187.0*ctx->mu*ctx->mu);
6419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->U,&y_ptr));
6429566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,ctx->U));
643c4762a1bSJed Brown 
644c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
6459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->U,&z_ptr));
6469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda[0],&y_ptr));
647c4762a1bSJed Brown   y_ptr[0] = 2.*(z_ptr[0]-ctx->ob[0]);
648c4762a1bSJed Brown   y_ptr[1] = 2.*(z_ptr[1]-ctx->ob[1]);
6499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda[0],&y_ptr));
6509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->U,&z_ptr));
6519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Mup[0],&y_ptr));
652c4762a1bSJed Brown   y_ptr[0] = 0.0;
6539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Mup[0],&y_ptr));
6549566063dSJacob Faibussowitsch   PetscCall(TSForwardGetSensitivities(ts,NULL,&tlmsen));
6559566063dSJacob Faibussowitsch   PetscCall(MatDenseGetColumn(tlmsen,0,&x_ptr));
6569566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda2[0],&y_ptr));
657c4762a1bSJed Brown   y_ptr[0] = 2.*x_ptr[0];
658c4762a1bSJed Brown   y_ptr[1] = 2.*x_ptr[1];
6599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0],&y_ptr));
6609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Mup2[0],&y_ptr));
661c4762a1bSJed Brown   y_ptr[0] = 0.0;
6629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Mup2[0],&y_ptr));
6639566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreColumn(tlmsen,&x_ptr));
6649566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts,1,ctx->Lambda,ctx->Mup));
665c4762a1bSJed Brown   if (ctx->implicitform) {
6669566063dSJacob Faibussowitsch     PetscCall(TSSetIHessianProduct(ts,ctx->Ihp1,IHessianProductUU,ctx->Ihp2,IHessianProductUP,ctx->Ihp3,IHessianProductPU,ctx->Ihp4,IHessianProductPP,ctx));
667c4762a1bSJed Brown   } else {
6689566063dSJacob Faibussowitsch     PetscCall(TSSetRHSHessianProduct(ts,ctx->Ihp1,RHSHessianProductUU,ctx->Ihp2,RHSHessianProductUP,ctx->Ihp3,RHSHessianProductPU,ctx->Ihp4,RHSHessianProductPP,ctx));
669c4762a1bSJed Brown   }
6709566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
671c4762a1bSJed Brown 
6729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda[0],&x_ptr));
6739566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->Lambda2[0],&y_ptr));
6749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->Mup2[0],&z_ptr));
675c4762a1bSJed Brown 
676c4762a1bSJed Brown   arr[0] = x_ptr[1]*dzdp2 + y_ptr[1]*dzdp2 + z_ptr[0];
677c4762a1bSJed Brown 
6789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0],&x_ptr));
6799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->Lambda2[0],&y_ptr));
6809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->Mup2[0],&z_ptr));
681c4762a1bSJed Brown 
682c4762a1bSJed Brown   /* Disable second-order adjoint mode */
6839566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
6849566063dSJacob Faibussowitsch   PetscCall(TSAdjointResetForward(ts));
685c4762a1bSJed Brown   PetscFunctionReturn(0);
686c4762a1bSJed Brown }
687c4762a1bSJed Brown 
688c4762a1bSJed Brown /*TEST
689c4762a1bSJed Brown     build:
690c4762a1bSJed Brown       requires: !complex !single
691c4762a1bSJed Brown     test:
692c4762a1bSJed 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
693c4762a1bSJed Brown       output_file: output/ex20opt_p_1.out
694c4762a1bSJed Brown 
695c4762a1bSJed Brown     test:
696c4762a1bSJed Brown       suffix: 2
697c4762a1bSJed 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
698c4762a1bSJed Brown       output_file: output/ex20opt_p_2.out
699c4762a1bSJed Brown 
700c4762a1bSJed Brown     test:
701c4762a1bSJed Brown       suffix: 3
702c4762a1bSJed Brown       args:  -ts_type cn -ts_adapt_type none -mu 100 -ts_dt 0.01 -viewer_binary_skip_info -tao_monitor -tao_view
703c4762a1bSJed Brown       output_file: output/ex20opt_p_3.out
704c4762a1bSJed Brown 
705c4762a1bSJed Brown     test:
706c4762a1bSJed Brown       suffix: 4
707c4762a1bSJed 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
708c4762a1bSJed Brown       output_file: output/ex20opt_p_4.out
709c4762a1bSJed Brown 
710c4762a1bSJed Brown TEST*/
711