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