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