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