xref: /petsc/src/ts/tutorials/ex20opt_ic.c (revision 0e3d61c972ee8b0cd7b6ee2ab64f8543b0740577)
1c4762a1bSJed Brown static char help[] = "Solves a ODE-constrained optimization problem -- finding the optimal initial conditions for the van der Pol equation.\n";
2c4762a1bSJed Brown 
3*0e3d61c9SBarry Smith /*
4c4762a1bSJed Brown   Concepts: TS^time-dependent nonlinear problems
5c4762a1bSJed Brown   Concepts: TS^van der Pol equation DAE equivalent
6ee300463SSatish Balay   Concepts: TS^Optimization using adjoint sensitivity analysis
7c4762a1bSJed Brown   Processors: 1
8c4762a1bSJed Brown */
9*0e3d61c9SBarry Smith /*
10c4762a1bSJed Brown   Notes:
11c4762a1bSJed Brown   This code demonstrates how to solve an ODE-constrained optimization problem with TAO, TSAdjoint and TS.
12c4762a1bSJed Brown   The nonlinear problem is written in an ODE equivalent form.
13c4762a1bSJed Brown   The objective is to minimize the difference between observation and model prediction by finding optimal values for initial conditions.
14c4762a1bSJed Brown   The gradient is computed with the discrete adjoint of an implicit method or an explicit method, see ex20adj.c for details.
15c4762a1bSJed Brown */
16c4762a1bSJed Brown 
17c4762a1bSJed Brown #include <petsctao.h>
18c4762a1bSJed Brown #include <petscts.h>
19c4762a1bSJed Brown 
20c4762a1bSJed Brown typedef struct _n_User *User;
21c4762a1bSJed Brown struct _n_User {
22c4762a1bSJed Brown   TS        ts;
23c4762a1bSJed Brown   PetscReal mu;
24c4762a1bSJed Brown   PetscReal next_output;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* Sensitivity analysis support */
27c4762a1bSJed Brown   PetscInt  steps;
28c4762a1bSJed Brown   PetscReal ftime;
29c4762a1bSJed Brown   Mat       A;                       /* Jacobian matrix for ODE */
30c4762a1bSJed Brown   Mat       Jacp;                    /* JacobianP matrix for ODE*/
31c4762a1bSJed Brown   Mat       H;                       /* Hessian matrix for optimization */
32c4762a1bSJed Brown   Vec       U,Lambda[1],Mup[1];      /* first-order adjoint variables */
33c4762a1bSJed Brown   Vec       Lambda2[2];              /* second-order adjoint variables */
34c4762a1bSJed Brown   Vec       Ihp1[1];                 /* working space for Hessian evaluations */
35c4762a1bSJed Brown   Vec       Dir;                     /* direction vector */
36c4762a1bSJed Brown   PetscReal ob[2];                   /* observation used by the cost function */
37c4762a1bSJed Brown   PetscBool implicitform;            /* implicit ODE? */
38c4762a1bSJed Brown };
39c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec,PetscScalar[],User);
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE  -------------------- */
42c4762a1bSJed Brown 
43c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
44c4762a1bSJed Brown {
45c4762a1bSJed Brown   PetscErrorCode    ierr;
46c4762a1bSJed Brown   User              user = (User)ctx;
47c4762a1bSJed Brown   PetscScalar       *f;
48c4762a1bSJed Brown   const PetscScalar *u;
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   PetscFunctionBeginUser;
51c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
52c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
53c4762a1bSJed Brown   f[0] = u[1];
54c4762a1bSJed Brown   f[1] = user->mu*((1.-u[0]*u[0])*u[1]-u[0]);
55c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
56c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
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   PetscErrorCode    ierr;
63c4762a1bSJed Brown   User              user = (User)ctx;
64c4762a1bSJed Brown   PetscReal         mu   = user->mu;
65c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
66c4762a1bSJed Brown   PetscScalar       J[2][2];
67c4762a1bSJed Brown   const PetscScalar *u;
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   PetscFunctionBeginUser;
70c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
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]);
75c4762a1bSJed Brown   ierr    = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78c4762a1bSJed Brown   if (A != B) {
79c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
81c4762a1bSJed Brown   }
82c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
83c4762a1bSJed Brown   PetscFunctionReturn(0);
84c4762a1bSJed Brown }
85c4762a1bSJed Brown 
86c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
87c4762a1bSJed Brown {
88c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
89c4762a1bSJed Brown   PetscScalar       *vhv;
90c4762a1bSJed Brown   PetscScalar       dJdU[2][2][2]={{{0}}};
91c4762a1bSJed Brown   PetscInt          i,j,k;
92c4762a1bSJed Brown   User              user = (User)ctx;
93c4762a1bSJed Brown   PetscErrorCode    ierr;
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   PetscFunctionBeginUser;
96c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
97c4762a1bSJed Brown   ierr = VecGetArrayRead(Vl[0],&vl);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = VecGetArrayRead(Vr,&vr);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = VecGetArray(VHV[0],&vhv);CHKERRQ(ierr);
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   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
111c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Vl[0],&vl);CHKERRQ(ierr);
112c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Vr,&vr);CHKERRQ(ierr);
113c4762a1bSJed Brown   ierr = VecRestoreArray(VHV[0],&vhv);CHKERRQ(ierr);
114c4762a1bSJed Brown   PetscFunctionReturn(0);
115c4762a1bSJed Brown }
116c4762a1bSJed Brown 
117c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE  -------------------- */
118c4762a1bSJed Brown 
119c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
120c4762a1bSJed Brown {
121c4762a1bSJed Brown   PetscErrorCode    ierr;
122c4762a1bSJed Brown   User              user = (User)ctx;
123c4762a1bSJed Brown   const PetscScalar *u,*udot;
124c4762a1bSJed Brown   PetscScalar       *f;
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   PetscFunctionBeginUser;
127c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
128c4762a1bSJed Brown   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
129c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
130c4762a1bSJed Brown   f[0] = udot[0] - u[1];
131c4762a1bSJed Brown   f[1] = udot[1] - user->mu*((1.0-u[0]*u[0])*u[1] - u[0]) ;
132c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
133c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
134c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
135c4762a1bSJed Brown   PetscFunctionReturn(0);
136c4762a1bSJed Brown }
137c4762a1bSJed Brown 
138c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
139c4762a1bSJed Brown {
140c4762a1bSJed Brown   PetscErrorCode    ierr;
141c4762a1bSJed Brown   User              user = (User)ctx;
142c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
143c4762a1bSJed Brown   PetscScalar       J[2][2];
144c4762a1bSJed Brown   const PetscScalar *u;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   PetscFunctionBeginUser;
147c4762a1bSJed Brown   ierr    = VecGetArrayRead(U,&u);CHKERRQ(ierr);
148c4762a1bSJed Brown   J[0][0] = a;     J[0][1] =  -1.0;
149c4762a1bSJed 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]);
150c4762a1bSJed Brown   ierr    = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
151c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
152c4762a1bSJed Brown   ierr    = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
153c4762a1bSJed Brown   ierr    = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154c4762a1bSJed Brown   if (A != B) {
155c4762a1bSJed Brown     ierr  = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
156c4762a1bSJed Brown     ierr  = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157c4762a1bSJed Brown   }
158c4762a1bSJed Brown   PetscFunctionReturn(0);
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
162c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx)
163c4762a1bSJed Brown {
164c4762a1bSJed Brown   PetscErrorCode    ierr;
165c4762a1bSJed Brown   const PetscScalar *u;
166c4762a1bSJed Brown   PetscReal         tfinal, dt;
167c4762a1bSJed Brown   User              user = (User)ctx;
168c4762a1bSJed Brown   Vec               interpolatedU;
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   PetscFunctionBeginUser;
171c4762a1bSJed Brown   ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
172c4762a1bSJed Brown   ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr);
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
175c4762a1bSJed Brown     ierr = VecDuplicate(U,&interpolatedU);CHKERRQ(ierr);
176c4762a1bSJed Brown     ierr = TSInterpolate(ts,user->next_output,interpolatedU);CHKERRQ(ierr);
177c4762a1bSJed Brown     ierr = VecGetArrayRead(interpolatedU,&u);CHKERRQ(ierr);
178c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"[%g] %D TS %g (dt = %g) X %g %g\n",
179c4762a1bSJed Brown                        (double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(u[0]),
180c4762a1bSJed Brown                        (double)PetscRealPart(u[1]));CHKERRQ(ierr);
181c4762a1bSJed Brown     ierr = VecRestoreArrayRead(interpolatedU,&u);CHKERRQ(ierr);
182c4762a1bSJed Brown     ierr = VecDestroy(&interpolatedU);CHKERRQ(ierr);
183c4762a1bSJed Brown     user->next_output += 0.1;
184c4762a1bSJed Brown   }
185c4762a1bSJed Brown   PetscFunctionReturn(0);
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown static PetscErrorCode IHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
189c4762a1bSJed Brown {
190c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
191c4762a1bSJed Brown   PetscScalar       *vhv;
192c4762a1bSJed Brown   PetscScalar       dJdU[2][2][2]={{{0}}};
193c4762a1bSJed Brown   PetscInt          i,j,k;
194c4762a1bSJed Brown   User              user = (User)ctx;
195c4762a1bSJed Brown   PetscErrorCode    ierr;
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   PetscFunctionBeginUser;
198c4762a1bSJed Brown   ierr          = VecGetArrayRead(U,&u);CHKERRQ(ierr);
199c4762a1bSJed Brown   ierr          = VecGetArrayRead(Vl[0],&vl);CHKERRQ(ierr);
200c4762a1bSJed Brown   ierr          = VecGetArrayRead(Vr,&vr);CHKERRQ(ierr);
201c4762a1bSJed Brown   ierr          = VecGetArray(VHV[0],&vhv);CHKERRQ(ierr);
202c4762a1bSJed Brown   dJdU[1][0][0] = 2.*user->mu*u[1];
203c4762a1bSJed Brown   dJdU[1][1][0] = 2.*user->mu*u[0];
204c4762a1bSJed Brown   dJdU[1][0][1] = 2.*user->mu*u[0];
205c4762a1bSJed Brown   for (j=0;j<2;j++) {
206c4762a1bSJed Brown     vhv[j] = 0;
207c4762a1bSJed Brown     for (k=0;k<2;k++)
208c4762a1bSJed Brown       for (i=0;i<2;i++)
209c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
210c4762a1bSJed Brown   }
211c4762a1bSJed Brown   ierr          = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
212c4762a1bSJed Brown   ierr          = VecRestoreArrayRead(Vl[0],&vl);CHKERRQ(ierr);
213c4762a1bSJed Brown   ierr          = VecRestoreArrayRead(Vr,&vr);CHKERRQ(ierr);
214c4762a1bSJed Brown   ierr          = VecRestoreArray(VHV[0],&vhv);CHKERRQ(ierr);
215c4762a1bSJed Brown   PetscFunctionReturn(0);
216c4762a1bSJed Brown }
217c4762a1bSJed Brown 
218c4762a1bSJed Brown /* ------------------ User-defined routines for TAO -------------------------- */
219c4762a1bSJed Brown 
220c4762a1bSJed Brown static PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx)
221c4762a1bSJed Brown {
222c4762a1bSJed Brown   User              user_ptr = (User)ctx;
223c4762a1bSJed Brown   TS                ts = user_ptr->ts;
224c4762a1bSJed Brown   const PetscScalar *x_ptr;
225c4762a1bSJed Brown   PetscScalar       *y_ptr;
226c4762a1bSJed Brown   PetscErrorCode ierr;
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   PetscFunctionBeginUser;
229c4762a1bSJed Brown   ierr = VecCopy(IC,user_ptr->U);CHKERRQ(ierr); /* set up the initial condition */
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   ierr = TSSetTime(ts,0.0);CHKERRQ(ierr);
232c4762a1bSJed Brown   ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr);
233c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,0.001);CHKERRQ(ierr); /* can be overwritten by command line options */
234c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = TSSolve(ts,user_ptr->U);CHKERRQ(ierr);
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   ierr     = VecGetArrayRead(user_ptr->U,&x_ptr);CHKERRQ(ierr);
238c4762a1bSJed Brown   ierr     = VecGetArray(user_ptr->Lambda[0],&y_ptr);CHKERRQ(ierr);
239c4762a1bSJed Brown   *f       = (x_ptr[0]-user_ptr->ob[0])*(x_ptr[0]-user_ptr->ob[0])+(x_ptr[1]-user_ptr->ob[1])*(x_ptr[1]-user_ptr->ob[1]);
240c4762a1bSJed Brown   y_ptr[0] = 2.*(x_ptr[0]-user_ptr->ob[0]);
241c4762a1bSJed Brown   y_ptr[1] = 2.*(x_ptr[1]-user_ptr->ob[1]);
242c4762a1bSJed Brown   ierr     = VecRestoreArray(user_ptr->Lambda[0],&y_ptr);CHKERRQ(ierr);
243c4762a1bSJed Brown   ierr     = VecRestoreArrayRead(user_ptr->U,&x_ptr);CHKERRQ(ierr);
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   ierr = TSSetCostGradients(ts,1,user_ptr->Lambda,NULL);CHKERRQ(ierr);
246c4762a1bSJed Brown   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
247c4762a1bSJed Brown   ierr = VecCopy(user_ptr->Lambda[0],G);CHKERRQ(ierr);
248c4762a1bSJed Brown   PetscFunctionReturn(0);
249c4762a1bSJed Brown }
250c4762a1bSJed Brown 
251c4762a1bSJed Brown static PetscErrorCode FormHessian(Tao tao,Vec U,Mat H,Mat Hpre,void *ctx)
252c4762a1bSJed Brown {
253c4762a1bSJed Brown   User           user_ptr = (User)ctx;
254c4762a1bSJed Brown   PetscScalar    harr[2];
255c4762a1bSJed Brown   PetscScalar    *x_ptr;
256c4762a1bSJed Brown   const PetscInt rows[2] = {0,1};
257c4762a1bSJed Brown   PetscInt       col;
258c4762a1bSJed Brown   PetscErrorCode ierr;
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   PetscFunctionBeginUser;
261c4762a1bSJed Brown   ierr     = VecCopy(U,user_ptr->U);CHKERRQ(ierr);
262c4762a1bSJed Brown   ierr     = VecGetArray(user_ptr->Dir,&x_ptr);CHKERRQ(ierr);
263c4762a1bSJed Brown   x_ptr[0] = 1.;
264c4762a1bSJed Brown   x_ptr[1] = 0.;
265c4762a1bSJed Brown   ierr     = VecRestoreArray(user_ptr->Dir,&x_ptr);CHKERRQ(ierr);
266c4762a1bSJed Brown   ierr     = Adjoint2(user_ptr->U,harr,user_ptr);CHKERRQ(ierr);
267c4762a1bSJed Brown   col      = 0;
268c4762a1bSJed Brown   ierr     = MatSetValues(H,2,rows,1,&col,harr,INSERT_VALUES);CHKERRQ(ierr);
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   ierr     = VecCopy(U,user_ptr->U);CHKERRQ(ierr);
271c4762a1bSJed Brown   ierr     = VecGetArray(user_ptr->Dir,&x_ptr);CHKERRQ(ierr);
272c4762a1bSJed Brown   x_ptr[0] = 0.;
273c4762a1bSJed Brown   x_ptr[1] = 1.;
274c4762a1bSJed Brown   ierr     = VecRestoreArray(user_ptr->Dir,&x_ptr);CHKERRQ(ierr);
275c4762a1bSJed Brown   ierr     = Adjoint2(user_ptr->U,harr,user_ptr);CHKERRQ(ierr);
276c4762a1bSJed Brown   col      = 1;
277c4762a1bSJed Brown   ierr     = MatSetValues(H,2,rows,1,&col,harr,INSERT_VALUES);CHKERRQ(ierr);
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   ierr   = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
280c4762a1bSJed Brown   ierr   = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
281c4762a1bSJed Brown   if (H != Hpre) {
282c4762a1bSJed Brown     ierr = MatAssemblyBegin(Hpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
283c4762a1bSJed Brown     ierr = MatAssemblyEnd(Hpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
284c4762a1bSJed Brown   }
285c4762a1bSJed Brown   PetscFunctionReturn(0);
286c4762a1bSJed Brown }
287c4762a1bSJed Brown 
288c4762a1bSJed Brown static PetscErrorCode MatrixFreeHessian(Tao tao,Vec U,Mat H,Mat Hpre,void *ctx)
289c4762a1bSJed Brown {
290c4762a1bSJed Brown   User           user_ptr = (User)ctx;
291c4762a1bSJed Brown   PetscErrorCode ierr;
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   PetscFunctionBeginUser;
294c4762a1bSJed Brown   ierr = VecCopy(U,user_ptr->U);CHKERRQ(ierr);
295c4762a1bSJed Brown   PetscFunctionReturn(0);
296c4762a1bSJed Brown }
297c4762a1bSJed Brown 
298c4762a1bSJed Brown /* ------------ Routines calculating second-order derivatives -------------- */
299c4762a1bSJed Brown 
300*0e3d61c9SBarry Smith /*
301c4762a1bSJed Brown   Compute the Hessian-vector product for the cost function using Second-order adjoint
302c4762a1bSJed Brown */
303c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec U,PetscScalar arr[],User ctx)
304c4762a1bSJed Brown {
305c4762a1bSJed Brown   TS             ts = ctx->ts;
306c4762a1bSJed Brown   PetscScalar    *x_ptr,*y_ptr;
307c4762a1bSJed Brown   Mat            tlmsen;
308c4762a1bSJed Brown   PetscErrorCode ierr;
309c4762a1bSJed Brown 
310c4762a1bSJed Brown   PetscFunctionBeginUser;
311c4762a1bSJed Brown   ierr = TSAdjointReset(ts);CHKERRQ(ierr);
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   ierr = TSSetTime(ts,0.0);CHKERRQ(ierr);
314c4762a1bSJed Brown   ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr);
315c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,0.001);CHKERRQ(ierr);
316c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
317c4762a1bSJed Brown   ierr = TSSetCostHessianProducts(ts,1,ctx->Lambda2,NULL,ctx->Dir);CHKERRQ(ierr);
318c4762a1bSJed Brown   ierr = TSAdjointSetForward(ts,NULL);CHKERRQ(ierr);
319c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
320c4762a1bSJed Brown 
321c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
322c4762a1bSJed Brown   ierr = VecGetArray(U,&x_ptr);CHKERRQ(ierr);
323c4762a1bSJed Brown   ierr = VecGetArray(ctx->Lambda[0],&y_ptr);CHKERRQ(ierr);
324c4762a1bSJed Brown   y_ptr[0] = 2.*(x_ptr[0]-ctx->ob[0]);
325c4762a1bSJed Brown   y_ptr[1] = 2.*(x_ptr[1]-ctx->ob[1]);
326c4762a1bSJed Brown   ierr = VecRestoreArray(ctx->Lambda[0],&y_ptr);CHKERRQ(ierr);
327c4762a1bSJed Brown   ierr = VecRestoreArray(U,&x_ptr);CHKERRQ(ierr);
328c4762a1bSJed Brown   ierr = TSForwardGetSensitivities(ts,NULL,&tlmsen);CHKERRQ(ierr);
329c4762a1bSJed Brown   ierr = MatDenseGetColumn(tlmsen,0,&x_ptr);CHKERRQ(ierr);
330c4762a1bSJed Brown   ierr = VecGetArray(ctx->Lambda2[0],&y_ptr);CHKERRQ(ierr);
331c4762a1bSJed Brown   y_ptr[0] = 2.*x_ptr[0];
332c4762a1bSJed Brown   y_ptr[1] = 2.*x_ptr[1];
333c4762a1bSJed Brown   ierr = VecRestoreArray(ctx->Lambda2[0],&y_ptr);CHKERRQ(ierr);
334c4762a1bSJed Brown   ierr = MatDenseRestoreColumn(tlmsen,&x_ptr);CHKERRQ(ierr);
335c4762a1bSJed Brown 
336c4762a1bSJed Brown   ierr = TSSetCostGradients(ts,1,ctx->Lambda,NULL);CHKERRQ(ierr);
337c4762a1bSJed Brown   if (ctx->implicitform) {
338c4762a1bSJed Brown     ierr = TSSetIHessianProduct(ts,ctx->Ihp1,IHessianProductUU,NULL,NULL,NULL,NULL,NULL,NULL,ctx);CHKERRQ(ierr);
339c4762a1bSJed Brown   } else {
340c4762a1bSJed Brown     ierr = TSSetRHSHessianProduct(ts,ctx->Ihp1,RHSHessianProductUU,NULL,NULL,NULL,NULL,NULL,NULL,ctx);CHKERRQ(ierr);
341c4762a1bSJed Brown   }
342c4762a1bSJed Brown   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
343c4762a1bSJed Brown 
344c4762a1bSJed Brown   ierr   = VecGetArray(ctx->Lambda2[0],&x_ptr);CHKERRQ(ierr);
345c4762a1bSJed Brown   arr[0] = x_ptr[0];
346c4762a1bSJed Brown   arr[1] = x_ptr[1];
347c4762a1bSJed Brown   ierr   = VecRestoreArray(ctx->Lambda2[0],&x_ptr);CHKERRQ(ierr);
348c4762a1bSJed Brown 
349c4762a1bSJed Brown   ierr = TSAdjointReset(ts);CHKERRQ(ierr);
350c4762a1bSJed Brown   ierr = TSAdjointResetForward(ts);CHKERRQ(ierr);
351c4762a1bSJed Brown   PetscFunctionReturn(0);
352c4762a1bSJed Brown }
353c4762a1bSJed Brown 
354c4762a1bSJed Brown PetscErrorCode FiniteDiff(Vec U,PetscScalar arr[],User ctx)
355c4762a1bSJed Brown {
356c4762a1bSJed Brown   Vec               Up,G,Gp;
357c4762a1bSJed Brown   const PetscScalar eps = PetscRealConstant(1e-7);
358c4762a1bSJed Brown   PetscScalar       *u;
359c4762a1bSJed Brown   Tao               tao = NULL;
360c4762a1bSJed Brown   PetscReal         f;
361c4762a1bSJed Brown   PetscErrorCode    ierr;
362c4762a1bSJed Brown 
363c4762a1bSJed Brown   PetscFunctionBeginUser;
364c4762a1bSJed Brown   ierr = VecDuplicate(U,&Up);CHKERRQ(ierr);
365c4762a1bSJed Brown   ierr = VecDuplicate(U,&G);CHKERRQ(ierr);
366c4762a1bSJed Brown   ierr = VecDuplicate(U,&Gp);CHKERRQ(ierr);
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   ierr = FormFunctionGradient(tao,U,&f,G,ctx);CHKERRQ(ierr);
369c4762a1bSJed Brown 
370c4762a1bSJed Brown   ierr = VecCopy(U,Up);CHKERRQ(ierr);
371c4762a1bSJed Brown   ierr = VecGetArray(Up,&u);CHKERRQ(ierr);
372c4762a1bSJed Brown   u[0] += eps;
373c4762a1bSJed Brown   ierr = VecRestoreArray(Up,&u);CHKERRQ(ierr);
374c4762a1bSJed Brown   ierr = FormFunctionGradient(tao,Up,&f,Gp,ctx);CHKERRQ(ierr);
375c4762a1bSJed Brown   ierr = VecAXPY(Gp,-1,G);CHKERRQ(ierr);
376c4762a1bSJed Brown   ierr = VecScale(Gp,1./eps);CHKERRQ(ierr);
377c4762a1bSJed Brown   ierr = VecGetArray(Gp,&u);CHKERRQ(ierr);
378c4762a1bSJed Brown   arr[0] = u[0];
379c4762a1bSJed Brown   arr[1] = u[1];
380c4762a1bSJed Brown   ierr  = VecRestoreArray(Gp,&u);CHKERRQ(ierr);
381c4762a1bSJed Brown 
382c4762a1bSJed Brown   ierr = VecCopy(U,Up);CHKERRQ(ierr);
383c4762a1bSJed Brown   ierr = VecGetArray(Up,&u);CHKERRQ(ierr);
384c4762a1bSJed Brown   u[1] += eps;
385c4762a1bSJed Brown   ierr = VecRestoreArray(Up,&u);CHKERRQ(ierr);
386c4762a1bSJed Brown   ierr = FormFunctionGradient(tao,Up,&f,Gp,ctx);CHKERRQ(ierr);
387c4762a1bSJed Brown   ierr = VecAXPY(Gp,-1,G);CHKERRQ(ierr);
388c4762a1bSJed Brown   ierr = VecScale(Gp,1./eps);CHKERRQ(ierr);
389c4762a1bSJed Brown   ierr = VecGetArray(Gp,&u);CHKERRQ(ierr);
390c4762a1bSJed Brown   arr[2] = u[0];
391c4762a1bSJed Brown   arr[3] = u[1];
392c4762a1bSJed Brown   ierr  = VecRestoreArray(Gp,&u);CHKERRQ(ierr);
393c4762a1bSJed Brown 
394c4762a1bSJed Brown   ierr = VecDestroy(&G);CHKERRQ(ierr);
395c4762a1bSJed Brown   ierr = VecDestroy(&Gp);CHKERRQ(ierr);
396c4762a1bSJed Brown   ierr = VecDestroy(&Up);CHKERRQ(ierr);
397c4762a1bSJed Brown   PetscFunctionReturn(0);
398c4762a1bSJed Brown }
399c4762a1bSJed Brown 
400c4762a1bSJed Brown static PetscErrorCode HessianProductMat(Mat mat,Vec svec,Vec y)
401c4762a1bSJed Brown {
402c4762a1bSJed Brown   User           user_ptr;
403c4762a1bSJed Brown   PetscScalar    *y_ptr;
404c4762a1bSJed Brown   PetscErrorCode ierr;
405c4762a1bSJed Brown 
406c4762a1bSJed Brown   PetscFunctionBeginUser;
407c4762a1bSJed Brown   ierr = MatShellGetContext(mat,(void*)&user_ptr);CHKERRQ(ierr);
408c4762a1bSJed Brown   ierr = VecCopy(svec,user_ptr->Dir);CHKERRQ(ierr);
409c4762a1bSJed Brown   ierr = VecGetArray(y,&y_ptr);CHKERRQ(ierr);
410c4762a1bSJed Brown   ierr = Adjoint2(user_ptr->U,y_ptr,user_ptr);CHKERRQ(ierr);
411c4762a1bSJed Brown   ierr = VecRestoreArray(y,&y_ptr);CHKERRQ(ierr);
412c4762a1bSJed Brown   PetscFunctionReturn(0);
413c4762a1bSJed Brown }
414c4762a1bSJed Brown 
415c4762a1bSJed Brown int main(int argc,char **argv)
416c4762a1bSJed Brown {
417c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE,mf = PETSC_TRUE;
418c4762a1bSJed Brown   PetscInt       mode = 0;
419c4762a1bSJed Brown   PetscMPIInt    size;
420c4762a1bSJed Brown   struct _n_User user;
421c4762a1bSJed Brown   Vec            x; /* working vector for TAO */
422c4762a1bSJed Brown   PetscScalar    *x_ptr,arr[4];
423c4762a1bSJed Brown   PetscScalar    ic1 = 2.2,ic2 = -0.7; /* initial guess for TAO */
424c4762a1bSJed Brown   Tao            tao;
425c4762a1bSJed Brown   KSP            ksp;
426c4762a1bSJed Brown   PC             pc;
427c4762a1bSJed Brown   PetscErrorCode ierr;
428c4762a1bSJed Brown 
429c4762a1bSJed Brown   /* Initialize program */
430c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
431ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
432c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
433c4762a1bSJed Brown 
434c4762a1bSJed Brown   /* Set runtime options */
435c4762a1bSJed Brown   user.next_output  = 0.0;
436c4762a1bSJed Brown   user.mu           = 1.0e3;
437c4762a1bSJed Brown   user.steps        = 0;
438c4762a1bSJed Brown   user.ftime        = 0.5;
439c4762a1bSJed Brown   user.implicitform = PETSC_TRUE;
440c4762a1bSJed Brown 
441c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr);
442c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr);
443c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-mode",&mode,NULL);CHKERRQ(ierr);
444c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-ic1",&ic1,NULL);CHKERRQ(ierr);
445c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-ic2",&ic2,NULL);CHKERRQ(ierr);
446c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-my_tao_mf",&mf,NULL);CHKERRQ(ierr); /* matrix-free hessian for optimization */
447c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-implicitform",&user.implicitform,NULL);CHKERRQ(ierr);
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
450c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&user.A);CHKERRQ(ierr);
451c4762a1bSJed Brown   ierr = MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
452c4762a1bSJed Brown   ierr = MatSetFromOptions(user.A);CHKERRQ(ierr);
453c4762a1bSJed Brown   ierr = MatSetUp(user.A);CHKERRQ(ierr);
454c4762a1bSJed Brown   ierr = MatCreateVecs(user.A,&user.U,NULL);CHKERRQ(ierr);
455c4762a1bSJed Brown   ierr = MatCreateVecs(user.A,&user.Dir,NULL);CHKERRQ(ierr);
456c4762a1bSJed Brown   ierr = MatCreateVecs(user.A,&user.Lambda[0],NULL);CHKERRQ(ierr);
457c4762a1bSJed Brown   ierr = MatCreateVecs(user.A,&user.Lambda2[0],NULL);CHKERRQ(ierr);
458c4762a1bSJed Brown   ierr = MatCreateVecs(user.A,&user.Ihp1[0],NULL);CHKERRQ(ierr);
459c4762a1bSJed Brown 
460c4762a1bSJed Brown   /* Create timestepping solver context */
461c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&user.ts);CHKERRQ(ierr);
462c4762a1bSJed Brown   ierr = TSSetEquationType(user.ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
463c4762a1bSJed Brown   if (user.implicitform) {
464c4762a1bSJed Brown     ierr = TSSetIFunction(user.ts,NULL,IFunction,&user);CHKERRQ(ierr);
465c4762a1bSJed Brown     ierr = TSSetIJacobian(user.ts,user.A,user.A,IJacobian,&user);CHKERRQ(ierr);
466c4762a1bSJed Brown     ierr = TSSetType(user.ts,TSCN);CHKERRQ(ierr);
467c4762a1bSJed Brown   } else {
468c4762a1bSJed Brown     ierr = TSSetRHSFunction(user.ts,NULL,RHSFunction,&user);CHKERRQ(ierr);
469c4762a1bSJed Brown     ierr = TSSetRHSJacobian(user.ts,user.A,user.A,RHSJacobian,&user);CHKERRQ(ierr);
470c4762a1bSJed Brown     ierr = TSSetType(user.ts,TSRK);CHKERRQ(ierr);
471c4762a1bSJed Brown   }
472c4762a1bSJed Brown   ierr = TSSetMaxTime(user.ts,user.ftime);CHKERRQ(ierr);
473c4762a1bSJed Brown   ierr = TSSetExactFinalTime(user.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
474c4762a1bSJed Brown 
475c4762a1bSJed Brown   if (monitor) {
476c4762a1bSJed Brown     ierr = TSMonitorSet(user.ts,Monitor,&user,NULL);CHKERRQ(ierr);
477c4762a1bSJed Brown   }
478c4762a1bSJed Brown 
479c4762a1bSJed Brown   /* Set ODE initial conditions */
480c4762a1bSJed Brown   ierr     = VecGetArray(user.U,&x_ptr);CHKERRQ(ierr);
481c4762a1bSJed Brown   x_ptr[0] = 2.0;
482c4762a1bSJed Brown   x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu);
483c4762a1bSJed Brown   ierr     = VecRestoreArray(user.U,&x_ptr);CHKERRQ(ierr);
484c4762a1bSJed Brown 
485c4762a1bSJed Brown   /* Set runtime options */
486c4762a1bSJed Brown   ierr = TSSetFromOptions(user.ts);CHKERRQ(ierr);
487c4762a1bSJed Brown 
488c4762a1bSJed Brown   /* Obtain the observation */
489c4762a1bSJed Brown   ierr       = TSSolve(user.ts,user.U);CHKERRQ(ierr);
490c4762a1bSJed Brown   ierr       = VecGetArray(user.U,&x_ptr);CHKERRQ(ierr);
491c4762a1bSJed Brown   user.ob[0] = x_ptr[0];
492c4762a1bSJed Brown   user.ob[1] = x_ptr[1];
493c4762a1bSJed Brown   ierr       = VecRestoreArray(user.U,&x_ptr);CHKERRQ(ierr);
494c4762a1bSJed Brown 
495c4762a1bSJed Brown   ierr     = VecDuplicate(user.U,&x);CHKERRQ(ierr);
496c4762a1bSJed Brown   ierr     = VecGetArray(x,&x_ptr);CHKERRQ(ierr);
497c4762a1bSJed Brown   x_ptr[0] = ic1;
498c4762a1bSJed Brown   x_ptr[1] = ic2;
499c4762a1bSJed Brown   ierr     = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr);
500c4762a1bSJed Brown 
501c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used */
502c4762a1bSJed Brown   ierr = TSSetSaveTrajectory(user.ts);CHKERRQ(ierr);
503c4762a1bSJed Brown 
504c4762a1bSJed Brown   /* Compare finite difference and second-order adjoint. */
505c4762a1bSJed Brown   switch (mode) {
506c4762a1bSJed Brown     case 2 :
507c4762a1bSJed Brown       ierr = FiniteDiff(x,arr,&user);CHKERRQ(ierr);
508c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD,"\n Finite difference approximation of the Hessian\n");CHKERRQ(ierr);
509c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD,"%g %g\n%g %g\n",(double)arr[0],(double)arr[1],(double)arr[2],(double)arr[3]);CHKERRQ(ierr);
510c4762a1bSJed Brown       break;
511c4762a1bSJed Brown     case 1 : /* Compute the Hessian column by column */
512c4762a1bSJed Brown       ierr     = VecCopy(x,user.U);CHKERRQ(ierr);
513c4762a1bSJed Brown       ierr     = VecGetArray(user.Dir,&x_ptr);CHKERRQ(ierr);
514c4762a1bSJed Brown       x_ptr[0] = 1.;
515c4762a1bSJed Brown       x_ptr[1] = 0.;
516c4762a1bSJed Brown       ierr     = VecRestoreArray(user.Dir,&x_ptr);CHKERRQ(ierr);
517c4762a1bSJed Brown       ierr     = Adjoint2(user.U,arr,&user);CHKERRQ(ierr);
518c4762a1bSJed Brown       ierr     = PetscPrintf(PETSC_COMM_WORLD,"\nFirst column of the Hessian\n");CHKERRQ(ierr);
519c4762a1bSJed Brown       ierr     = PetscPrintf(PETSC_COMM_WORLD,"%g\n%g\n",(double)arr[0],(double)arr[1]);CHKERRQ(ierr);
520c4762a1bSJed Brown       ierr     = VecCopy(x,user.U);CHKERRQ(ierr);
521c4762a1bSJed Brown       ierr     = VecGetArray(user.Dir,&x_ptr);CHKERRQ(ierr);
522c4762a1bSJed Brown       x_ptr[0] = 0.;
523c4762a1bSJed Brown       x_ptr[1] = 1.;
524c4762a1bSJed Brown       ierr     = VecRestoreArray(user.Dir,&x_ptr);CHKERRQ(ierr);
525c4762a1bSJed Brown       ierr     = Adjoint2(user.U,arr,&user);CHKERRQ(ierr);
526c4762a1bSJed Brown       ierr     = PetscPrintf(PETSC_COMM_WORLD,"\nSecond column of the Hessian\n");CHKERRQ(ierr);
527c4762a1bSJed Brown       ierr     = PetscPrintf(PETSC_COMM_WORLD,"%g\n%g\n",(double)arr[0],(double)arr[1]);CHKERRQ(ierr);
528c4762a1bSJed Brown       break;
529c4762a1bSJed Brown     case 0 :
530c4762a1bSJed Brown       /* Create optimization context and set up */
531c4762a1bSJed Brown       ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
532c4762a1bSJed Brown       ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr);
533c4762a1bSJed Brown       ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&user);CHKERRQ(ierr);
534c4762a1bSJed Brown 
535c4762a1bSJed Brown       if (mf) {
536c4762a1bSJed Brown         ierr = MatCreateShell(PETSC_COMM_SELF,2,2,2,2,(void*)&user,&user.H);CHKERRQ(ierr);
537c4762a1bSJed Brown         ierr = MatShellSetOperation(user.H,MATOP_MULT,(void(*)(void))HessianProductMat);CHKERRQ(ierr);
538c4762a1bSJed Brown         ierr = MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
539c4762a1bSJed Brown         ierr = TaoSetHessianRoutine(tao,user.H,user.H,MatrixFreeHessian,(void *)&user);CHKERRQ(ierr);
540c4762a1bSJed Brown       } else { /* Create Hessian matrix */
541c4762a1bSJed Brown         ierr = MatCreate(PETSC_COMM_WORLD,&user.H);CHKERRQ(ierr);
542c4762a1bSJed Brown         ierr = MatSetSizes(user.H,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
543c4762a1bSJed Brown         ierr = MatSetUp(user.H);CHKERRQ(ierr);
544c4762a1bSJed Brown         ierr = TaoSetHessianRoutine(tao,user.H,user.H,FormHessian,(void *)&user);CHKERRQ(ierr);
545c4762a1bSJed Brown       }
546c4762a1bSJed Brown 
547c4762a1bSJed Brown       /* Not use any preconditioner */
548c4762a1bSJed Brown       ierr   = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
549c4762a1bSJed Brown       if (ksp) {
550c4762a1bSJed Brown         ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
551c4762a1bSJed Brown         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
552c4762a1bSJed Brown       }
553c4762a1bSJed Brown 
554c4762a1bSJed Brown       /* Set initial solution guess */
555c4762a1bSJed Brown       ierr = TaoSetInitialVector(tao,x);CHKERRQ(ierr);
556c4762a1bSJed Brown       ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
557c4762a1bSJed Brown       ierr = TaoSolve(tao);CHKERRQ(ierr);
558c4762a1bSJed Brown       ierr = TaoDestroy(&tao);CHKERRQ(ierr);
559c4762a1bSJed Brown       ierr = MatDestroy(&user.H);CHKERRQ(ierr);
560c4762a1bSJed Brown       break;
561c4762a1bSJed Brown     default:
562c4762a1bSJed Brown       break;
563c4762a1bSJed Brown   }
564c4762a1bSJed Brown 
565c4762a1bSJed Brown   /* Free work space.  All PETSc objects should be destroyed when they are no longer needed. */
566c4762a1bSJed Brown   ierr = MatDestroy(&user.A);CHKERRQ(ierr);
567c4762a1bSJed Brown   ierr = VecDestroy(&user.U);CHKERRQ(ierr);
568c4762a1bSJed Brown   ierr = VecDestroy(&user.Lambda[0]);CHKERRQ(ierr);
569c4762a1bSJed Brown   ierr = VecDestroy(&user.Lambda2[0]);CHKERRQ(ierr);
570c4762a1bSJed Brown   ierr = VecDestroy(&user.Ihp1[0]);CHKERRQ(ierr);
571c4762a1bSJed Brown   ierr = VecDestroy(&user.Dir);CHKERRQ(ierr);
572c4762a1bSJed Brown   ierr = TSDestroy(&user.ts);CHKERRQ(ierr);
573c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
574c4762a1bSJed Brown   ierr = PetscFinalize();
575c4762a1bSJed Brown   return(ierr);
576c4762a1bSJed Brown }
577c4762a1bSJed Brown 
578c4762a1bSJed Brown /*TEST
579c4762a1bSJed Brown     build:
580c4762a1bSJed Brown       requires: !complex !single
581c4762a1bSJed Brown 
582c4762a1bSJed Brown     test:
583c4762a1bSJed Brown       args:  -ts_type cn -viewer_binary_skip_info -tao_monitor -tao_view -mu 1000 -ts_dt 0.03125
584c4762a1bSJed Brown       output_file: output/ex20opt_ic_1.out
585c4762a1bSJed Brown 
586c4762a1bSJed Brown     test:
587c4762a1bSJed Brown       suffix: 2
588c4762a1bSJed Brown       args:  -ts_type beuler -viewer_binary_skip_info -tao_monitor -tao_view -mu 100 -ts_dt 0.01 -tao_type bntr -tao_bnk_pc_type none
589c4762a1bSJed Brown       output_file: output/ex20opt_ic_2.out
590c4762a1bSJed Brown 
591c4762a1bSJed Brown     test:
592c4762a1bSJed Brown       suffix: 3
593c4762a1bSJed Brown       args:  -ts_type cn -viewer_binary_skip_info -tao_monitor -tao_view -mu 100 -ts_dt 0.01 -tao_type bntr -tao_bnk_pc_type none
594c4762a1bSJed Brown       output_file: output/ex20opt_ic_3.out
595c4762a1bSJed Brown 
596c4762a1bSJed Brown     test:
597c4762a1bSJed Brown       suffix: 4
598c4762a1bSJed Brown       args: -implicitform 0 -ts_dt 0.01 -ts_max_steps 2 -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view -mode 1 -my_tao_mf
599c4762a1bSJed Brown TEST*/
600