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