xref: /petsc/src/ts/tutorials/ex20opt_p.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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    Concepts: TS^time-dependent nonlinear problems
7c4762a1bSJed Brown    Concepts: TS^van der Pol equation DAE equivalent
8ee300463SSatish Balay    Concepts: TS^Optimization using adjoint sensitivity analysis
9c4762a1bSJed Brown    Processors: 1
10c4762a1bSJed Brown */
11c4762a1bSJed Brown /* ------------------------------------------------------------------------
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   Notes:
14c4762a1bSJed Brown   This code demonstrates how to solve a DAE-constrained optimization problem with TAO, TSAdjoint and TS.
15c4762a1bSJed Brown   The nonlinear problem is written in a DAE equivalent form.
16c4762a1bSJed Brown   The objective is to minimize the difference between observation and model prediction by finding an optimal value for parameter \mu.
17c4762a1bSJed Brown   The gradient is computed with the discrete adjoint of an implicit theta method, see ex20adj.c for details.
18c4762a1bSJed Brown   ------------------------------------------------------------------------- */
19c4762a1bSJed Brown #include <petsctao.h>
20c4762a1bSJed Brown #include <petscts.h>
21c4762a1bSJed Brown 
22c4762a1bSJed Brown typedef struct _n_User *User;
23c4762a1bSJed Brown struct _n_User {
24c4762a1bSJed Brown   TS        ts;
25c4762a1bSJed Brown   PetscReal mu;
26c4762a1bSJed Brown   PetscReal next_output;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Sensitivity analysis support */
29c4762a1bSJed Brown   PetscReal ftime;
30c4762a1bSJed Brown   Mat       A;                       /* Jacobian matrix */
31c4762a1bSJed Brown   Mat       Jacp;                    /* JacobianP matrix */
32c4762a1bSJed Brown   Mat       H;                       /* Hessian matrix for optimization */
33c4762a1bSJed Brown   Vec       U,Lambda[1],Mup[1];      /* adjoint variables */
34c4762a1bSJed Brown   Vec       Lambda2[1],Mup2[1];      /* second-order adjoint variables */
35c4762a1bSJed Brown   Vec       Ihp1[1];                 /* working space for Hessian evaluations */
36c4762a1bSJed Brown   Vec       Ihp2[1];                 /* working space for Hessian evaluations */
37c4762a1bSJed Brown   Vec       Ihp3[1];                 /* working space for Hessian evaluations */
38c4762a1bSJed Brown   Vec       Ihp4[1];                 /* working space for Hessian evaluations */
39c4762a1bSJed Brown   Vec       Dir;                     /* direction vector */
40c4762a1bSJed Brown   PetscReal ob[2];                   /* observation used by the cost function */
41c4762a1bSJed Brown   PetscBool implicitform;            /* implicit ODE? */
42c4762a1bSJed Brown };
43c4762a1bSJed Brown 
44c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
45c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
46c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec,PetscScalar[],User);
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE  -------------------- */
49c4762a1bSJed Brown 
50c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
51c4762a1bSJed Brown {
52c4762a1bSJed Brown   User              user = (User)ctx;
53c4762a1bSJed Brown   PetscScalar       *f;
54c4762a1bSJed Brown   const PetscScalar *u;
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   PetscFunctionBeginUser;
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
59c4762a1bSJed Brown   f[0] = u[1];
60c4762a1bSJed Brown   f[1] = user->mu*((1.-u[0]*u[0])*u[1]-u[0]);
61*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
62*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
63c4762a1bSJed Brown   PetscFunctionReturn(0);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
67c4762a1bSJed Brown {
68c4762a1bSJed Brown   User              user = (User)ctx;
69c4762a1bSJed Brown   PetscReal         mu   = user->mu;
70c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
71c4762a1bSJed Brown   PetscScalar       J[2][2];
72c4762a1bSJed Brown   const PetscScalar *u;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   PetscFunctionBeginUser;
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   J[0][0] = 0;
78c4762a1bSJed Brown   J[1][0] = -mu*(2.0*u[1]*u[0]+1.);
79c4762a1bSJed Brown   J[0][1] = 1.0;
80c4762a1bSJed Brown   J[1][1] = mu*(1.0-u[0]*u[0]);
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
84c4762a1bSJed Brown   if (B && A != B) {
85*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
86*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
87c4762a1bSJed Brown   }
88c4762a1bSJed Brown 
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
90c4762a1bSJed Brown   PetscFunctionReturn(0);
91c4762a1bSJed Brown }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
94c4762a1bSJed Brown {
95c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
96c4762a1bSJed Brown   PetscScalar       *vhv;
97c4762a1bSJed Brown   PetscScalar       dJdU[2][2][2]={{{0}}};
98c4762a1bSJed Brown   PetscInt          i,j,k;
99c4762a1bSJed Brown   User              user = (User)ctx;
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   PetscFunctionBeginUser;
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
103*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   dJdU[1][0][0] = -2.*user->mu*u[1];
108c4762a1bSJed Brown   dJdU[1][1][0] = -2.*user->mu*u[0];
109c4762a1bSJed Brown   dJdU[1][0][1] = -2.*user->mu*u[0];
110c4762a1bSJed Brown   for (j=0; j<2; j++) {
111c4762a1bSJed Brown     vhv[j] = 0;
112c4762a1bSJed Brown     for (k=0; k<2; k++)
113c4762a1bSJed Brown       for (i=0; i<2; i++)
114c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
115c4762a1bSJed Brown   }
116c4762a1bSJed Brown 
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
121c4762a1bSJed Brown   PetscFunctionReturn(0);
122c4762a1bSJed Brown }
123c4762a1bSJed Brown 
124c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
125c4762a1bSJed Brown {
126c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
127c4762a1bSJed Brown   PetscScalar       *vhv;
128c4762a1bSJed Brown   PetscScalar       dJdP[2][2][1]={{{0}}};
129c4762a1bSJed Brown   PetscInt          i,j,k;
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   PetscFunctionBeginUser;
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   dJdP[1][0][0] = -(1.+2.*u[0]*u[1]);
138c4762a1bSJed Brown   dJdP[1][1][0] = 1.-u[0]*u[0];
139c4762a1bSJed Brown   for (j=0; j<2; j++) {
140c4762a1bSJed Brown     vhv[j] = 0;
141c4762a1bSJed Brown     for (k=0; k<1; k++)
142c4762a1bSJed Brown       for (i=0; i<2; i++)
143c4762a1bSJed Brown         vhv[j] += vl[i]*dJdP[i][j][k]*vr[k];
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
150c4762a1bSJed Brown   PetscFunctionReturn(0);
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
154c4762a1bSJed Brown {
155c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
156c4762a1bSJed Brown   PetscScalar       *vhv;
157c4762a1bSJed Brown   PetscScalar       dJdU[2][1][2]={{{0}}};
158c4762a1bSJed Brown   PetscInt          i,j,k;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   PetscFunctionBeginUser;
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
164*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   dJdU[1][0][0] = -1.-2.*u[1]*u[0];
167c4762a1bSJed Brown   dJdU[1][0][1] = 1.-u[0]*u[0];
168c4762a1bSJed Brown   for (j=0; j<1; j++) {
169c4762a1bSJed Brown     vhv[j] = 0;
170c4762a1bSJed Brown     for (k=0; k<2; k++)
171c4762a1bSJed Brown       for (i=0; i<2; i++)
172c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
173c4762a1bSJed Brown   }
174c4762a1bSJed Brown 
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
176*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
177*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
178*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
179c4762a1bSJed Brown   PetscFunctionReturn(0);
180c4762a1bSJed Brown }
181c4762a1bSJed Brown 
182c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
183c4762a1bSJed Brown {
184c4762a1bSJed Brown   PetscFunctionBeginUser;
185c4762a1bSJed Brown   PetscFunctionReturn(0);
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE  -------------------- */
189c4762a1bSJed Brown 
190c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
191c4762a1bSJed Brown {
192c4762a1bSJed Brown   User              user = (User)ctx;
193c4762a1bSJed Brown   PetscScalar       *f;
194c4762a1bSJed Brown   const PetscScalar *u,*udot;
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   PetscFunctionBeginUser;
197*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
198*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Udot,&udot));
199*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   f[0] = udot[0] - u[1];
202c4762a1bSJed Brown   f[1] = udot[1] - user->mu*((1.0-u[0]*u[0])*u[1] - u[0]) ;
203c4762a1bSJed Brown 
204*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
205*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Udot,&udot));
206*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
207c4762a1bSJed Brown   PetscFunctionReturn(0);
208c4762a1bSJed Brown }
209c4762a1bSJed Brown 
210c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
211c4762a1bSJed Brown {
212c4762a1bSJed Brown   User              user     = (User)ctx;
213c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
214c4762a1bSJed Brown   PetscScalar       J[2][2];
215c4762a1bSJed Brown   const PetscScalar *u;
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   PetscFunctionBeginUser;
218*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   J[0][0] = a;     J[0][1] =  -1.0;
221c4762a1bSJed 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]);
222*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
223*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
224*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
225*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
226c4762a1bSJed Brown   if (A != B) {
227*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
228*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
229c4762a1bSJed Brown   }
230c4762a1bSJed Brown 
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
232c4762a1bSJed Brown   PetscFunctionReturn(0);
233c4762a1bSJed Brown }
234c4762a1bSJed Brown 
235c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec U,Mat A,void *ctx)
236c4762a1bSJed Brown {
237c4762a1bSJed Brown   PetscInt          row[] = {0,1},col[]={0};
238c4762a1bSJed Brown   PetscScalar       J[2][1];
239c4762a1bSJed Brown   const PetscScalar *u;
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   PetscFunctionBeginUser;
242*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   J[0][0] = 0;
245c4762a1bSJed Brown   J[1][0] = (1.-u[0]*u[0])*u[1]-u[0];
246*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES));
247*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
248*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
249*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
250c4762a1bSJed Brown 
251*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
252c4762a1bSJed Brown   PetscFunctionReturn(0);
253c4762a1bSJed Brown }
254c4762a1bSJed Brown 
255c4762a1bSJed Brown static PetscErrorCode IHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
256c4762a1bSJed Brown {
257c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
258c4762a1bSJed Brown   PetscScalar       *vhv;
259c4762a1bSJed Brown   PetscScalar       dJdU[2][2][2]={{{0}}};
260c4762a1bSJed Brown   PetscInt          i,j,k;
261c4762a1bSJed Brown   User              user = (User)ctx;
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   PetscFunctionBeginUser;
264*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
265*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
266*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
267*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   dJdU[1][0][0] = 2.*user->mu*u[1];
270c4762a1bSJed Brown   dJdU[1][1][0] = 2.*user->mu*u[0];
271c4762a1bSJed Brown   dJdU[1][0][1] = 2.*user->mu*u[0];
272c4762a1bSJed Brown   for (j=0; j<2; j++) {
273c4762a1bSJed Brown     vhv[j] = 0;
274c4762a1bSJed Brown     for (k=0; k<2; k++)
275c4762a1bSJed Brown       for (i=0; i<2; i++)
276c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
277c4762a1bSJed Brown   }
278c4762a1bSJed Brown 
279*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
280*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
281*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
282*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
283c4762a1bSJed Brown   PetscFunctionReturn(0);
284c4762a1bSJed Brown }
285c4762a1bSJed Brown 
286c4762a1bSJed Brown static PetscErrorCode IHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
287c4762a1bSJed Brown {
288c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
289c4762a1bSJed Brown   PetscScalar       *vhv;
290c4762a1bSJed Brown   PetscScalar       dJdP[2][2][1]={{{0}}};
291c4762a1bSJed Brown   PetscInt          i,j,k;
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   PetscFunctionBeginUser;
294*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
295*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
296*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
297*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
298c4762a1bSJed Brown 
299c4762a1bSJed Brown   dJdP[1][0][0] = 1.+2.*u[0]*u[1];
300c4762a1bSJed Brown   dJdP[1][1][0] = u[0]*u[0]-1.;
301c4762a1bSJed Brown   for (j=0; j<2; j++) {
302c4762a1bSJed Brown     vhv[j] = 0;
303c4762a1bSJed Brown     for (k=0; k<1; k++)
304c4762a1bSJed Brown       for (i=0; i<2; i++)
305c4762a1bSJed Brown         vhv[j] += vl[i]*dJdP[i][j][k]*vr[k];
306c4762a1bSJed Brown   }
307c4762a1bSJed Brown 
308*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
309*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
310*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
311*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
312c4762a1bSJed Brown   PetscFunctionReturn(0);
313c4762a1bSJed Brown }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown static PetscErrorCode IHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
316c4762a1bSJed Brown {
317c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
318c4762a1bSJed Brown   PetscScalar       *vhv;
319c4762a1bSJed Brown   PetscScalar       dJdU[2][1][2]={{{0}}};
320c4762a1bSJed Brown   PetscInt          i,j,k;
321c4762a1bSJed Brown 
322c4762a1bSJed Brown   PetscFunctionBeginUser;
323*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
324*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
325*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
326*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   dJdU[1][0][0] = 1.+2.*u[1]*u[0];
329c4762a1bSJed Brown   dJdU[1][0][1] = u[0]*u[0]-1.;
330c4762a1bSJed Brown   for (j=0; j<1; j++) {
331c4762a1bSJed Brown     vhv[j] = 0;
332c4762a1bSJed Brown     for (k=0; k<2; k++)
333c4762a1bSJed Brown       for (i=0; i<2; i++)
334c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
335c4762a1bSJed Brown   }
336c4762a1bSJed Brown 
337*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
338*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
339*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
340*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
341c4762a1bSJed Brown   PetscFunctionReturn(0);
342c4762a1bSJed Brown }
343c4762a1bSJed Brown 
344c4762a1bSJed Brown static PetscErrorCode IHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
345c4762a1bSJed Brown {
346c4762a1bSJed Brown   PetscFunctionBeginUser;
347c4762a1bSJed Brown   PetscFunctionReturn(0);
348c4762a1bSJed Brown }
349c4762a1bSJed Brown 
350c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
351c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
352c4762a1bSJed Brown {
353c4762a1bSJed Brown   PetscErrorCode    ierr;
354c4762a1bSJed Brown   const PetscScalar *x;
355c4762a1bSJed Brown   PetscReal         tfinal, dt;
356c4762a1bSJed Brown   User              user = (User)ctx;
357c4762a1bSJed Brown   Vec               interpolatedX;
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   PetscFunctionBeginUser;
360*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTimeStep(ts,&dt));
361*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetMaxTime(ts,&tfinal));
362c4762a1bSJed Brown 
363c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
364*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(X,&interpolatedX));
365*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSInterpolate(ts,user->next_output,interpolatedX));
366*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(interpolatedX,&x));
367c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"[%g] %D TS %g (dt = %g) X %g %g\n",
368c4762a1bSJed Brown                        (double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(x[0]),
369c4762a1bSJed Brown                        (double)PetscRealPart(x[1]));CHKERRQ(ierr);
370*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(interpolatedX,&x));
371*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&interpolatedX));
372c4762a1bSJed Brown     user->next_output += PetscRealConstant(0.1);
373c4762a1bSJed Brown   }
374c4762a1bSJed Brown   PetscFunctionReturn(0);
375c4762a1bSJed Brown }
376c4762a1bSJed Brown 
377c4762a1bSJed Brown int main(int argc,char **argv)
378c4762a1bSJed Brown {
379c4762a1bSJed Brown   Vec                P;
380c4762a1bSJed Brown   PetscBool          monitor = PETSC_FALSE;
381c4762a1bSJed Brown   PetscScalar        *x_ptr;
382c4762a1bSJed Brown   const PetscScalar  *y_ptr;
383c4762a1bSJed Brown   PetscMPIInt        size;
384c4762a1bSJed Brown   struct _n_User     user;
385c4762a1bSJed Brown   PetscErrorCode     ierr;
386c4762a1bSJed Brown   Tao                tao;
387c4762a1bSJed Brown   KSP                ksp;
388c4762a1bSJed Brown   PC                 pc;
389c4762a1bSJed Brown 
390c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
391c4762a1bSJed Brown      Initialize program
392c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
393c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
394*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
3953c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
396c4762a1bSJed Brown 
397c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
398*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
399*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOBQNLS));
400c4762a1bSJed Brown 
401c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
402c4762a1bSJed Brown     Set runtime options
403c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
404c4762a1bSJed Brown   user.next_output  = 0.0;
405c4762a1bSJed Brown   user.mu           = PetscRealConstant(1.0e3);
406c4762a1bSJed Brown   user.ftime        = PetscRealConstant(0.5);
407c4762a1bSJed Brown   user.implicitform = PETSC_TRUE;
408c4762a1bSJed Brown 
409*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
410*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
411*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-implicitform",&user.implicitform,NULL));
412c4762a1bSJed Brown 
413c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
414*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.A));
415*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2));
416*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user.A));
417*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user.A));
418*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.U,NULL));
419*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.Lambda[0],NULL));
420*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.Lambda2[0],NULL));
421*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.Ihp1[0],NULL));
422*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.Ihp2[0],NULL));
423c4762a1bSJed Brown 
424*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.Jacp));
425*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1));
426*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user.Jacp));
427*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user.Jacp));
428*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&user.Dir,NULL));
429*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&user.Mup[0],NULL));
430*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&user.Mup2[0],NULL));
431*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&user.Ihp3[0],NULL));
432*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&user.Ihp4[0],NULL));
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   /* Create timestepping solver context */
435*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&user.ts));
436*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetEquationType(user.ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
437c4762a1bSJed Brown   if (user.implicitform) {
438*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetIFunction(user.ts,NULL,IFunction,&user));
439*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetIJacobian(user.ts,user.A,user.A,IJacobian,&user));
440*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetType(user.ts,TSCN));
441c4762a1bSJed Brown   } else {
442*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSFunction(user.ts,NULL,RHSFunction,&user));
443*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSJacobian(user.ts,user.A,user.A,RHSJacobian,&user));
444*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetType(user.ts,TSRK));
445c4762a1bSJed Brown   }
446*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobianP(user.ts,user.Jacp,RHSJacobianP,&user));
447*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(user.ts,user.ftime));
448*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(user.ts,TS_EXACTFINALTIME_MATCHSTEP));
449c4762a1bSJed Brown 
450c4762a1bSJed Brown   if (monitor) {
451*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorSet(user.ts,Monitor,&user,NULL));
452c4762a1bSJed Brown   }
453c4762a1bSJed Brown 
454c4762a1bSJed Brown   /* Set ODE initial conditions */
455*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user.U,&x_ptr));
456c4762a1bSJed Brown   x_ptr[0] = 2.0;
457c4762a1bSJed Brown   x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu);
458*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user.U,&x_ptr));
459*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(user.ts,PetscRealConstant(0.001)));
460c4762a1bSJed Brown 
461c4762a1bSJed Brown   /* Set runtime options */
462*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(user.ts));
463c4762a1bSJed Brown 
464*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(user.ts,user.U));
465*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(user.U,&y_ptr));
466c4762a1bSJed Brown   user.ob[0] = y_ptr[0];
467c4762a1bSJed Brown   user.ob[1] = y_ptr[1];
468*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(user.U,&y_ptr));
469c4762a1bSJed Brown 
470c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used.
471c4762a1bSJed Brown      Skip checkpointing for the first TSSolve since no adjoint run follows it.
472c4762a1bSJed Brown    */
473*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSaveTrajectory(user.ts));
474c4762a1bSJed Brown 
475c4762a1bSJed Brown   /* Optimization starts */
476*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.H));
477*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user.H,PETSC_DECIDE,PETSC_DECIDE,1,1));
478*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user.H)); /* Hessian should be symmetric. Do we need to do MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE) ? */
479c4762a1bSJed Brown 
480c4762a1bSJed Brown   /* Set initial solution guess */
481*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&P,NULL));
482*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(P,&x_ptr));
483c4762a1bSJed Brown   x_ptr[0] = PetscRealConstant(1.2);
484*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(P,&x_ptr));
485*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,P));
486c4762a1bSJed Brown 
487c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
488*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user));
489*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user));
490c4762a1bSJed Brown 
491c4762a1bSJed Brown   /* Check for any TAO command line options */
492*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetKSP(tao,&ksp));
493c4762a1bSJed Brown   if (ksp) {
494*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp,&pc));
495*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetType(pc,PCNONE));
496c4762a1bSJed Brown   }
497*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
498c4762a1bSJed Brown 
499*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
500c4762a1bSJed Brown 
501*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(P,PETSC_VIEWER_STDOUT_WORLD));
502c4762a1bSJed Brown   /* Free TAO data structures */
503*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
504c4762a1bSJed Brown 
505c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
506c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
507c4762a1bSJed Brown      are no longer needed.
508c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
509*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.H));
510*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.A));
511*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.U));
512*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.Jacp));
513*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Lambda[0]));
514*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Mup[0]));
515*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Lambda2[0]));
516*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Mup2[0]));
517*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Ihp1[0]));
518*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Ihp2[0]));
519*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Ihp3[0]));
520*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Ihp4[0]));
521*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Dir));
522*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&user.ts));
523*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&P));
524c4762a1bSJed Brown   ierr = PetscFinalize();
525c4762a1bSJed Brown   return ierr;
526c4762a1bSJed Brown }
527c4762a1bSJed Brown 
528c4762a1bSJed Brown /* ------------------------------------------------------------------ */
529c4762a1bSJed Brown /*
530c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
531c4762a1bSJed Brown 
532c4762a1bSJed Brown    Input Parameters:
533c4762a1bSJed Brown    tao - the Tao context
534c4762a1bSJed Brown    X   - the input vector
535a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
536c4762a1bSJed Brown 
537c4762a1bSJed Brown    Output Parameters:
538c4762a1bSJed Brown    f   - the newly evaluated function
539c4762a1bSJed Brown    G   - the newly evaluated gradient
540c4762a1bSJed Brown */
541c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
542c4762a1bSJed Brown {
543c4762a1bSJed Brown   User              user_ptr = (User)ctx;
544c4762a1bSJed Brown   TS                ts = user_ptr->ts;
545c4762a1bSJed Brown   PetscScalar       *x_ptr,*g;
546c4762a1bSJed Brown   const PetscScalar *y_ptr;
547c4762a1bSJed Brown 
548c4762a1bSJed Brown   PetscFunctionBeginUser;
549*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,&y_ptr));
550c4762a1bSJed Brown   user_ptr->mu = y_ptr[0];
551*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,&y_ptr));
552c4762a1bSJed Brown 
553*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0.0));
554*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ts,0));
555*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,PetscRealConstant(0.001))); /* can be overwritten by command line options */
556*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
557*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user_ptr->U,&x_ptr));
558c4762a1bSJed Brown   x_ptr[0] = 2.0;
559c4762a1bSJed 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);
560*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user_ptr->U,&x_ptr));
561c4762a1bSJed Brown 
562*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,user_ptr->U));
563c4762a1bSJed Brown 
564*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(user_ptr->U,&y_ptr));
565c4762a1bSJed 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]);
566c4762a1bSJed Brown 
567c4762a1bSJed Brown   /*   Reset initial conditions for the adjoint integration */
568*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user_ptr->Lambda[0],&x_ptr));
569c4762a1bSJed Brown   x_ptr[0] = 2.*(y_ptr[0]-user_ptr->ob[0]);
570c4762a1bSJed Brown   x_ptr[1] = 2.*(y_ptr[1]-user_ptr->ob[1]);
571*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(user_ptr->U,&y_ptr));
572*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user_ptr->Lambda[0],&x_ptr));
573c4762a1bSJed Brown 
574*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user_ptr->Mup[0],&x_ptr));
575c4762a1bSJed Brown   x_ptr[0] = 0.0;
576*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user_ptr->Mup[0],&x_ptr));
577*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,1,user_ptr->Lambda,user_ptr->Mup));
578c4762a1bSJed Brown 
579*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSolve(ts));
580c4762a1bSJed Brown 
581*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user_ptr->Mup[0],&x_ptr));
582*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(user_ptr->Lambda[0],&y_ptr));
583*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(G,&g));
584c4762a1bSJed 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];
585*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user_ptr->Mup[0],&x_ptr));
586*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(user_ptr->Lambda[0],&y_ptr));
587*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(G,&g));
588c4762a1bSJed Brown   PetscFunctionReturn(0);
589c4762a1bSJed Brown }
590c4762a1bSJed Brown 
591c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec P,Mat H,Mat Hpre,void *ctx)
592c4762a1bSJed Brown {
593c4762a1bSJed Brown   User           user_ptr = (User)ctx;
594c4762a1bSJed Brown   PetscScalar    harr[1];
595c4762a1bSJed Brown   const PetscInt rows[1] = {0};
596c4762a1bSJed Brown   PetscInt       col = 0;
597c4762a1bSJed Brown 
598c4762a1bSJed Brown   PetscFunctionBeginUser;
599*5f80ce2aSJacob Faibussowitsch   CHKERRQ(Adjoint2(P,harr,user_ptr));
600*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(H,1,rows,1,&col,harr,INSERT_VALUES));
601c4762a1bSJed Brown 
602*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
603*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
604c4762a1bSJed Brown   if (H != Hpre) {
605*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(Hpre,MAT_FINAL_ASSEMBLY));
606*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(Hpre,MAT_FINAL_ASSEMBLY));
607c4762a1bSJed Brown   }
608c4762a1bSJed Brown   PetscFunctionReturn(0);
609c4762a1bSJed Brown }
610c4762a1bSJed Brown 
611c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec P,PetscScalar arr[],User ctx)
612c4762a1bSJed Brown {
613c4762a1bSJed Brown   TS                ts = ctx->ts;
614c4762a1bSJed Brown   const PetscScalar *z_ptr;
615c4762a1bSJed Brown   PetscScalar       *x_ptr,*y_ptr,dzdp,dzdp2;
616c4762a1bSJed Brown   Mat               tlmsen;
617c4762a1bSJed Brown 
618c4762a1bSJed Brown   PetscFunctionBeginUser;
619c4762a1bSJed Brown   /* Reset TSAdjoint so that AdjointSetUp will be called again */
620*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointReset(ts));
621c4762a1bSJed Brown 
622c4762a1bSJed Brown   /* The directional vector should be 1 since it is one-dimensional */
623*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Dir,&x_ptr));
624c4762a1bSJed Brown   x_ptr[0] = 1.;
625*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Dir,&x_ptr));
626c4762a1bSJed Brown 
627*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,&z_ptr));
628c4762a1bSJed Brown   ctx->mu = z_ptr[0];
629*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,&z_ptr));
630c4762a1bSJed Brown 
631c4762a1bSJed Brown   dzdp  = -10.0/(81.0*ctx->mu*ctx->mu) + 2.0*292.0/(2187.0*ctx->mu*ctx->mu*ctx->mu);
632c4762a1bSJed 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);
633c4762a1bSJed Brown 
634*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0.0));
635*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ts,0));
636*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,PetscRealConstant(0.001))); /* can be overwritten by command line options */
637*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
638*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostHessianProducts(ts,1,ctx->Lambda2,ctx->Mup2,ctx->Dir));
639c4762a1bSJed Brown 
640*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(ctx->Jacp));
641*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(ctx->Jacp,1,0,dzdp,INSERT_VALUES));
642*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(ctx->Jacp,MAT_FINAL_ASSEMBLY));
643*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(ctx->Jacp,MAT_FINAL_ASSEMBLY));
644c4762a1bSJed Brown 
645*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSetForward(ts,ctx->Jacp));
646*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->U,&y_ptr));
647c4762a1bSJed Brown   y_ptr[0] = 2.0;
648c4762a1bSJed Brown   y_ptr[1] = -2.0/3.0 + 10.0/(81.0*ctx->mu) - 292.0/(2187.0*ctx->mu*ctx->mu);
649*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->U,&y_ptr));
650*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,ctx->U));
651c4762a1bSJed Brown 
652c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
653*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(ctx->U,&z_ptr));
654*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Lambda[0],&y_ptr));
655c4762a1bSJed Brown   y_ptr[0] = 2.*(z_ptr[0]-ctx->ob[0]);
656c4762a1bSJed Brown   y_ptr[1] = 2.*(z_ptr[1]-ctx->ob[1]);
657*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Lambda[0],&y_ptr));
658*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(ctx->U,&z_ptr));
659*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Mup[0],&y_ptr));
660c4762a1bSJed Brown   y_ptr[0] = 0.0;
661*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Mup[0],&y_ptr));
662*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSForwardGetSensitivities(ts,NULL,&tlmsen));
663*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetColumn(tlmsen,0,&x_ptr));
664*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Lambda2[0],&y_ptr));
665c4762a1bSJed Brown   y_ptr[0] = 2.*x_ptr[0];
666c4762a1bSJed Brown   y_ptr[1] = 2.*x_ptr[1];
667*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Lambda2[0],&y_ptr));
668*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Mup2[0],&y_ptr));
669c4762a1bSJed Brown   y_ptr[0] = 0.0;
670*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Mup2[0],&y_ptr));
671*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreColumn(tlmsen,&x_ptr));
672*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,1,ctx->Lambda,ctx->Mup));
673c4762a1bSJed Brown   if (ctx->implicitform) {
674*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetIHessianProduct(ts,ctx->Ihp1,IHessianProductUU,ctx->Ihp2,IHessianProductUP,ctx->Ihp3,IHessianProductPU,ctx->Ihp4,IHessianProductPP,ctx));
675c4762a1bSJed Brown   } else {
676*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSHessianProduct(ts,ctx->Ihp1,RHSHessianProductUU,ctx->Ihp2,RHSHessianProductUP,ctx->Ihp3,RHSHessianProductPU,ctx->Ihp4,RHSHessianProductPP,ctx));
677c4762a1bSJed Brown   }
678*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSolve(ts));
679c4762a1bSJed Brown 
680*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Lambda[0],&x_ptr));
681*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Lambda2[0],&y_ptr));
682*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(ctx->Mup2[0],&z_ptr));
683c4762a1bSJed Brown 
684c4762a1bSJed Brown   arr[0] = x_ptr[1]*dzdp2 + y_ptr[1]*dzdp2 + z_ptr[0];
685c4762a1bSJed Brown 
686*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Lambda2[0],&x_ptr));
687*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Lambda2[0],&y_ptr));
688*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(ctx->Mup2[0],&z_ptr));
689c4762a1bSJed Brown 
690c4762a1bSJed Brown   /* Disable second-order adjoint mode */
691*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointReset(ts));
692*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointResetForward(ts));
693c4762a1bSJed Brown   PetscFunctionReturn(0);
694c4762a1bSJed Brown }
695c4762a1bSJed Brown 
696c4762a1bSJed Brown /*TEST
697c4762a1bSJed Brown     build:
698c4762a1bSJed Brown       requires: !complex !single
699c4762a1bSJed Brown     test:
700c4762a1bSJed 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
701c4762a1bSJed Brown       output_file: output/ex20opt_p_1.out
702c4762a1bSJed Brown 
703c4762a1bSJed Brown     test:
704c4762a1bSJed Brown       suffix: 2
705c4762a1bSJed 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
706c4762a1bSJed Brown       output_file: output/ex20opt_p_2.out
707c4762a1bSJed Brown 
708c4762a1bSJed Brown     test:
709c4762a1bSJed Brown       suffix: 3
710c4762a1bSJed Brown       args:  -ts_type cn -ts_adapt_type none -mu 100 -ts_dt 0.01 -viewer_binary_skip_info -tao_monitor -tao_view
711c4762a1bSJed Brown       output_file: output/ex20opt_p_3.out
712c4762a1bSJed Brown 
713c4762a1bSJed Brown     test:
714c4762a1bSJed Brown       suffix: 4
715c4762a1bSJed 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
716c4762a1bSJed Brown       output_file: output/ex20opt_p_4.out
717c4762a1bSJed Brown 
718c4762a1bSJed Brown TEST*/
719