xref: /petsc/src/ts/tutorials/ex20opt_p.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation.\n\
3c4762a1bSJed Brown Input parameters include:\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown    Concepts: TS^time-dependent nonlinear problems
7c4762a1bSJed Brown    Concepts: TS^van der Pol equation DAE equivalent
8ee300463SSatish Balay    Concepts: TS^Optimization using adjoint sensitivity analysis
9c4762a1bSJed Brown    Processors: 1
10c4762a1bSJed Brown */
11c4762a1bSJed Brown /* ------------------------------------------------------------------------
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   Notes:
14c4762a1bSJed Brown   This code demonstrates how to solve a DAE-constrained optimization problem with TAO, TSAdjoint and TS.
15c4762a1bSJed Brown   The nonlinear problem is written in a DAE equivalent form.
16c4762a1bSJed Brown   The objective is to minimize the difference between observation and model prediction by finding an optimal value for parameter \mu.
17c4762a1bSJed Brown   The gradient is computed with the discrete adjoint of an implicit theta method, see ex20adj.c for details.
18c4762a1bSJed Brown   ------------------------------------------------------------------------- */
19c4762a1bSJed Brown #include <petsctao.h>
20c4762a1bSJed Brown #include <petscts.h>
21c4762a1bSJed Brown 
22c4762a1bSJed Brown typedef struct _n_User *User;
23c4762a1bSJed Brown struct _n_User {
24c4762a1bSJed Brown   TS        ts;
25c4762a1bSJed Brown   PetscReal mu;
26c4762a1bSJed Brown   PetscReal next_output;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Sensitivity analysis support */
29c4762a1bSJed Brown   PetscReal ftime;
30c4762a1bSJed Brown   Mat       A;                       /* Jacobian matrix */
31c4762a1bSJed Brown   Mat       Jacp;                    /* JacobianP matrix */
32c4762a1bSJed Brown   Mat       H;                       /* Hessian matrix for optimization */
33c4762a1bSJed Brown   Vec       U,Lambda[1],Mup[1];      /* adjoint variables */
34c4762a1bSJed Brown   Vec       Lambda2[1],Mup2[1];      /* second-order adjoint variables */
35c4762a1bSJed Brown   Vec       Ihp1[1];                 /* working space for Hessian evaluations */
36c4762a1bSJed Brown   Vec       Ihp2[1];                 /* working space for Hessian evaluations */
37c4762a1bSJed Brown   Vec       Ihp3[1];                 /* working space for Hessian evaluations */
38c4762a1bSJed Brown   Vec       Ihp4[1];                 /* working space for Hessian evaluations */
39c4762a1bSJed Brown   Vec       Dir;                     /* direction vector */
40c4762a1bSJed Brown   PetscReal ob[2];                   /* observation used by the cost function */
41c4762a1bSJed Brown   PetscBool implicitform;            /* implicit ODE? */
42c4762a1bSJed Brown };
43c4762a1bSJed Brown 
44c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
45c4762a1bSJed Brown PetscErrorCode FormHessian(Tao,Vec,Mat,Mat,void*);
46c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec,PetscScalar[],User);
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE  -------------------- */
49c4762a1bSJed Brown 
50c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
51c4762a1bSJed Brown {
52c4762a1bSJed Brown   User              user = (User)ctx;
53c4762a1bSJed Brown   PetscScalar       *f;
54c4762a1bSJed Brown   const PetscScalar *u;
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   PetscFunctionBeginUser;
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
59c4762a1bSJed Brown   f[0] = u[1];
60c4762a1bSJed Brown   f[1] = user->mu*((1.-u[0]*u[0])*u[1]-u[0]);
615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
63c4762a1bSJed Brown   PetscFunctionReturn(0);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
67c4762a1bSJed Brown {
68c4762a1bSJed Brown   User              user = (User)ctx;
69c4762a1bSJed Brown   PetscReal         mu   = user->mu;
70c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
71c4762a1bSJed Brown   PetscScalar       J[2][2];
72c4762a1bSJed Brown   const PetscScalar *u;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   PetscFunctionBeginUser;
755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   J[0][0] = 0;
78c4762a1bSJed Brown   J[1][0] = -mu*(2.0*u[1]*u[0]+1.);
79c4762a1bSJed Brown   J[0][1] = 1.0;
80c4762a1bSJed Brown   J[1][1] = mu*(1.0-u[0]*u[0]);
815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
84c4762a1bSJed Brown   if (B && A != B) {
855f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
865f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
87c4762a1bSJed Brown   }
88c4762a1bSJed Brown 
895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
90c4762a1bSJed Brown   PetscFunctionReturn(0);
91c4762a1bSJed Brown }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
94c4762a1bSJed Brown {
95c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
96c4762a1bSJed Brown   PetscScalar       *vhv;
97c4762a1bSJed Brown   PetscScalar       dJdU[2][2][2]={{{0}}};
98c4762a1bSJed Brown   PetscInt          i,j,k;
99c4762a1bSJed Brown   User              user = (User)ctx;
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   PetscFunctionBeginUser;
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   dJdU[1][0][0] = -2.*user->mu*u[1];
108c4762a1bSJed Brown   dJdU[1][1][0] = -2.*user->mu*u[0];
109c4762a1bSJed Brown   dJdU[1][0][1] = -2.*user->mu*u[0];
110c4762a1bSJed Brown   for (j=0; j<2; j++) {
111c4762a1bSJed Brown     vhv[j] = 0;
112c4762a1bSJed Brown     for (k=0; k<2; k++)
113c4762a1bSJed Brown       for (i=0; i<2; i++)
114c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
115c4762a1bSJed Brown   }
116c4762a1bSJed Brown 
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
121c4762a1bSJed Brown   PetscFunctionReturn(0);
122c4762a1bSJed Brown }
123c4762a1bSJed Brown 
124c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
125c4762a1bSJed Brown {
126c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
127c4762a1bSJed Brown   PetscScalar       *vhv;
128c4762a1bSJed Brown   PetscScalar       dJdP[2][2][1]={{{0}}};
129c4762a1bSJed Brown   PetscInt          i,j,k;
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   PetscFunctionBeginUser;
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   dJdP[1][0][0] = -(1.+2.*u[0]*u[1]);
138c4762a1bSJed Brown   dJdP[1][1][0] = 1.-u[0]*u[0];
139c4762a1bSJed Brown   for (j=0; j<2; j++) {
140c4762a1bSJed Brown     vhv[j] = 0;
141c4762a1bSJed Brown     for (k=0; k<1; k++)
142c4762a1bSJed Brown       for (i=0; i<2; i++)
143c4762a1bSJed Brown         vhv[j] += vl[i]*dJdP[i][j][k]*vr[k];
144c4762a1bSJed Brown   }
145c4762a1bSJed Brown 
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
150c4762a1bSJed Brown   PetscFunctionReturn(0);
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
154c4762a1bSJed Brown {
155c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
156c4762a1bSJed Brown   PetscScalar       *vhv;
157c4762a1bSJed Brown   PetscScalar       dJdU[2][1][2]={{{0}}};
158c4762a1bSJed Brown   PetscInt          i,j,k;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   PetscFunctionBeginUser;
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   dJdU[1][0][0] = -1.-2.*u[1]*u[0];
167c4762a1bSJed Brown   dJdU[1][0][1] = 1.-u[0]*u[0];
168c4762a1bSJed Brown   for (j=0; j<1; j++) {
169c4762a1bSJed Brown     vhv[j] = 0;
170c4762a1bSJed Brown     for (k=0; k<2; k++)
171c4762a1bSJed Brown       for (i=0; i<2; i++)
172c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
173c4762a1bSJed Brown   }
174c4762a1bSJed Brown 
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
179c4762a1bSJed Brown   PetscFunctionReturn(0);
180c4762a1bSJed Brown }
181c4762a1bSJed Brown 
182c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
183c4762a1bSJed Brown {
184c4762a1bSJed Brown   PetscFunctionBeginUser;
185c4762a1bSJed Brown   PetscFunctionReturn(0);
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE  -------------------- */
189c4762a1bSJed Brown 
190c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
191c4762a1bSJed Brown {
192c4762a1bSJed Brown   User              user = (User)ctx;
193c4762a1bSJed Brown   PetscScalar       *f;
194c4762a1bSJed Brown   const PetscScalar *u,*udot;
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   PetscFunctionBeginUser;
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Udot,&udot));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   f[0] = udot[0] - u[1];
202c4762a1bSJed Brown   f[1] = udot[1] - user->mu*((1.0-u[0]*u[0])*u[1] - u[0]) ;
203c4762a1bSJed Brown 
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Udot,&udot));
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
207c4762a1bSJed Brown   PetscFunctionReturn(0);
208c4762a1bSJed Brown }
209c4762a1bSJed Brown 
210c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
211c4762a1bSJed Brown {
212c4762a1bSJed Brown   User              user     = (User)ctx;
213c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
214c4762a1bSJed Brown   PetscScalar       J[2][2];
215c4762a1bSJed Brown   const PetscScalar *u;
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   PetscFunctionBeginUser;
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   J[0][0] = a;     J[0][1] =  -1.0;
221c4762a1bSJed 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]);
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
2245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
2255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
226c4762a1bSJed Brown   if (A != B) {
2275f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
2285f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
229c4762a1bSJed Brown   }
230c4762a1bSJed Brown 
2315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
232c4762a1bSJed Brown   PetscFunctionReturn(0);
233c4762a1bSJed Brown }
234c4762a1bSJed Brown 
235c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec U,Mat A,void *ctx)
236c4762a1bSJed Brown {
237c4762a1bSJed Brown   PetscInt          row[] = {0,1},col[]={0};
238c4762a1bSJed Brown   PetscScalar       J[2][1];
239c4762a1bSJed Brown   const PetscScalar *u;
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   PetscFunctionBeginUser;
2425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   J[0][0] = 0;
245c4762a1bSJed Brown   J[1][0] = (1.-u[0]*u[0])*u[1]-u[0];
2465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES));
2475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
250c4762a1bSJed Brown 
2515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
252c4762a1bSJed Brown   PetscFunctionReturn(0);
253c4762a1bSJed Brown }
254c4762a1bSJed Brown 
255c4762a1bSJed Brown static PetscErrorCode IHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
256c4762a1bSJed Brown {
257c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
258c4762a1bSJed Brown   PetscScalar       *vhv;
259c4762a1bSJed Brown   PetscScalar       dJdU[2][2][2]={{{0}}};
260c4762a1bSJed Brown   PetscInt          i,j,k;
261c4762a1bSJed Brown   User              user = (User)ctx;
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   PetscFunctionBeginUser;
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   dJdU[1][0][0] = 2.*user->mu*u[1];
270c4762a1bSJed Brown   dJdU[1][1][0] = 2.*user->mu*u[0];
271c4762a1bSJed Brown   dJdU[1][0][1] = 2.*user->mu*u[0];
272c4762a1bSJed Brown   for (j=0; j<2; j++) {
273c4762a1bSJed Brown     vhv[j] = 0;
274c4762a1bSJed Brown     for (k=0; k<2; k++)
275c4762a1bSJed Brown       for (i=0; i<2; i++)
276c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
277c4762a1bSJed Brown   }
278c4762a1bSJed Brown 
2795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
2805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
2815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
2825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
283c4762a1bSJed Brown   PetscFunctionReturn(0);
284c4762a1bSJed Brown }
285c4762a1bSJed Brown 
286c4762a1bSJed Brown static PetscErrorCode IHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
287c4762a1bSJed Brown {
288c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
289c4762a1bSJed Brown   PetscScalar       *vhv;
290c4762a1bSJed Brown   PetscScalar       dJdP[2][2][1]={{{0}}};
291c4762a1bSJed Brown   PetscInt          i,j,k;
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   PetscFunctionBeginUser;
2945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
2955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
2965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
2975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
298c4762a1bSJed Brown 
299c4762a1bSJed Brown   dJdP[1][0][0] = 1.+2.*u[0]*u[1];
300c4762a1bSJed Brown   dJdP[1][1][0] = u[0]*u[0]-1.;
301c4762a1bSJed Brown   for (j=0; j<2; j++) {
302c4762a1bSJed Brown     vhv[j] = 0;
303c4762a1bSJed Brown     for (k=0; k<1; k++)
304c4762a1bSJed Brown       for (i=0; i<2; i++)
305c4762a1bSJed Brown         vhv[j] += vl[i]*dJdP[i][j][k]*vr[k];
306c4762a1bSJed Brown   }
307c4762a1bSJed Brown 
3085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
3095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
3105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
3115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
312c4762a1bSJed Brown   PetscFunctionReturn(0);
313c4762a1bSJed Brown }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown static PetscErrorCode IHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
316c4762a1bSJed Brown {
317c4762a1bSJed Brown   const PetscScalar *vl,*vr,*u;
318c4762a1bSJed Brown   PetscScalar       *vhv;
319c4762a1bSJed Brown   PetscScalar       dJdU[2][1][2]={{{0}}};
320c4762a1bSJed Brown   PetscInt          i,j,k;
321c4762a1bSJed Brown 
322c4762a1bSJed Brown   PetscFunctionBeginUser;
3235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
3245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
3255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
3265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(VHV[0],&vhv));
327c4762a1bSJed Brown 
328c4762a1bSJed Brown   dJdU[1][0][0] = 1.+2.*u[1]*u[0];
329c4762a1bSJed Brown   dJdU[1][0][1] = u[0]*u[0]-1.;
330c4762a1bSJed Brown   for (j=0; j<1; j++) {
331c4762a1bSJed Brown     vhv[j] = 0;
332c4762a1bSJed Brown     for (k=0; k<2; k++)
333c4762a1bSJed Brown       for (i=0; i<2; i++)
334c4762a1bSJed Brown         vhv[j] += vl[i]*dJdU[i][j][k]*vr[k];
335c4762a1bSJed Brown   }
336c4762a1bSJed Brown 
3375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
3385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
3395f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
3405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(VHV[0],&vhv));
341c4762a1bSJed Brown   PetscFunctionReturn(0);
342c4762a1bSJed Brown }
343c4762a1bSJed Brown 
344c4762a1bSJed Brown static PetscErrorCode IHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
345c4762a1bSJed Brown {
346c4762a1bSJed Brown   PetscFunctionBeginUser;
347c4762a1bSJed Brown   PetscFunctionReturn(0);
348c4762a1bSJed Brown }
349c4762a1bSJed Brown 
350c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
351c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
352c4762a1bSJed Brown {
353c4762a1bSJed Brown   PetscErrorCode    ierr;
354c4762a1bSJed Brown   const PetscScalar *x;
355c4762a1bSJed Brown   PetscReal         tfinal, dt;
356c4762a1bSJed Brown   User              user = (User)ctx;
357c4762a1bSJed Brown   Vec               interpolatedX;
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   PetscFunctionBeginUser;
3605f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTimeStep(ts,&dt));
3615f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetMaxTime(ts,&tfinal));
362c4762a1bSJed Brown 
363c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
3645f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(X,&interpolatedX));
3655f80ce2aSJacob Faibussowitsch     CHKERRQ(TSInterpolate(ts,user->next_output,interpolatedX));
3665f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArrayRead(interpolatedX,&x));
367c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"[%g] %D TS %g (dt = %g) X %g %g\n",
368c4762a1bSJed Brown                        (double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(x[0]),
369c4762a1bSJed Brown                        (double)PetscRealPart(x[1]));CHKERRQ(ierr);
3705f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArrayRead(interpolatedX,&x));
3715f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&interpolatedX));
372c4762a1bSJed Brown     user->next_output += PetscRealConstant(0.1);
373c4762a1bSJed Brown   }
374c4762a1bSJed Brown   PetscFunctionReturn(0);
375c4762a1bSJed Brown }
376c4762a1bSJed Brown 
377c4762a1bSJed Brown int main(int argc,char **argv)
378c4762a1bSJed Brown {
379c4762a1bSJed Brown   Vec                P;
380c4762a1bSJed Brown   PetscBool          monitor = PETSC_FALSE;
381c4762a1bSJed Brown   PetscScalar        *x_ptr;
382c4762a1bSJed Brown   const PetscScalar  *y_ptr;
383c4762a1bSJed Brown   PetscMPIInt        size;
384c4762a1bSJed Brown   struct _n_User     user;
385c4762a1bSJed Brown   Tao                tao;
386c4762a1bSJed Brown   KSP                ksp;
387c4762a1bSJed Brown   PC                 pc;
388c4762a1bSJed Brown 
389c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
390c4762a1bSJed Brown      Initialize program
391c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
392*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
3935f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
3943c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
395c4762a1bSJed Brown 
396c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
3975f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
3985f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOBQNLS));
399c4762a1bSJed Brown 
400c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
401c4762a1bSJed Brown     Set runtime options
402c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
403c4762a1bSJed Brown   user.next_output  = 0.0;
404c4762a1bSJed Brown   user.mu           = PetscRealConstant(1.0e3);
405c4762a1bSJed Brown   user.ftime        = PetscRealConstant(0.5);
406c4762a1bSJed Brown   user.implicitform = PETSC_TRUE;
407c4762a1bSJed Brown 
4085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
4095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
4105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-implicitform",&user.implicitform,NULL));
411c4762a1bSJed Brown 
412c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
4135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.A));
4145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2));
4155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user.A));
4165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user.A));
4175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.U,NULL));
4185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.Lambda[0],NULL));
4195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.Lambda2[0],NULL));
4205f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.Ihp1[0],NULL));
4215f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.Ihp2[0],NULL));
422c4762a1bSJed Brown 
4235f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.Jacp));
4245f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1));
4255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user.Jacp));
4265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user.Jacp));
4275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&user.Dir,NULL));
4285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&user.Mup[0],NULL));
4295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&user.Mup2[0],NULL));
4305f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&user.Ihp3[0],NULL));
4315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&user.Ihp4[0],NULL));
432c4762a1bSJed Brown 
433c4762a1bSJed Brown   /* Create timestepping solver context */
4345f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&user.ts));
4355f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetEquationType(user.ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
436c4762a1bSJed Brown   if (user.implicitform) {
4375f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetIFunction(user.ts,NULL,IFunction,&user));
4385f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetIJacobian(user.ts,user.A,user.A,IJacobian,&user));
4395f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetType(user.ts,TSCN));
440c4762a1bSJed Brown   } else {
4415f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSFunction(user.ts,NULL,RHSFunction,&user));
4425f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSJacobian(user.ts,user.A,user.A,RHSJacobian,&user));
4435f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetType(user.ts,TSRK));
444c4762a1bSJed Brown   }
4455f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobianP(user.ts,user.Jacp,RHSJacobianP,&user));
4465f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(user.ts,user.ftime));
4475f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(user.ts,TS_EXACTFINALTIME_MATCHSTEP));
448c4762a1bSJed Brown 
449c4762a1bSJed Brown   if (monitor) {
4505f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorSet(user.ts,Monitor,&user,NULL));
451c4762a1bSJed Brown   }
452c4762a1bSJed Brown 
453c4762a1bSJed Brown   /* Set ODE initial conditions */
4545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user.U,&x_ptr));
455c4762a1bSJed Brown   x_ptr[0] = 2.0;
456c4762a1bSJed Brown   x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu);
4575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user.U,&x_ptr));
4585f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(user.ts,PetscRealConstant(0.001)));
459c4762a1bSJed Brown 
460c4762a1bSJed Brown   /* Set runtime options */
4615f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(user.ts));
462c4762a1bSJed Brown 
4635f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(user.ts,user.U));
4645f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(user.U,&y_ptr));
465c4762a1bSJed Brown   user.ob[0] = y_ptr[0];
466c4762a1bSJed Brown   user.ob[1] = y_ptr[1];
4675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(user.U,&y_ptr));
468c4762a1bSJed Brown 
469c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used.
470c4762a1bSJed Brown      Skip checkpointing for the first TSSolve since no adjoint run follows it.
471c4762a1bSJed Brown    */
4725f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSaveTrajectory(user.ts));
473c4762a1bSJed Brown 
474c4762a1bSJed Brown   /* Optimization starts */
4755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.H));
4765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user.H,PETSC_DECIDE,PETSC_DECIDE,1,1));
4775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user.H)); /* Hessian should be symmetric. Do we need to do MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE) ? */
478c4762a1bSJed Brown 
479c4762a1bSJed Brown   /* Set initial solution guess */
4805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&P,NULL));
4815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(P,&x_ptr));
482c4762a1bSJed Brown   x_ptr[0] = PetscRealConstant(1.2);
4835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(P,&x_ptr));
4845f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,P));
485c4762a1bSJed Brown 
486c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
4875f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user));
4885f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user));
489c4762a1bSJed Brown 
490c4762a1bSJed Brown   /* Check for any TAO command line options */
4915f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetKSP(tao,&ksp));
492c4762a1bSJed Brown   if (ksp) {
4935f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp,&pc));
4945f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetType(pc,PCNONE));
495c4762a1bSJed Brown   }
4965f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
497c4762a1bSJed Brown 
4985f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
499c4762a1bSJed Brown 
5005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(P,PETSC_VIEWER_STDOUT_WORLD));
501c4762a1bSJed Brown   /* Free TAO data structures */
5025f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
503c4762a1bSJed Brown 
504c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
505c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
506c4762a1bSJed Brown      are no longer needed.
507c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
5085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.H));
5095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.A));
5105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.U));
5115f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.Jacp));
5125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Lambda[0]));
5135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Mup[0]));
5145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Lambda2[0]));
5155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Mup2[0]));
5165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Ihp1[0]));
5175f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Ihp2[0]));
5185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Ihp3[0]));
5195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Ihp4[0]));
5205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.Dir));
5215f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&user.ts));
5225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&P));
523*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
524*b122ec5aSJacob Faibussowitsch   return 0;
525c4762a1bSJed Brown }
526c4762a1bSJed Brown 
527c4762a1bSJed Brown /* ------------------------------------------------------------------ */
528c4762a1bSJed Brown /*
529c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
530c4762a1bSJed Brown 
531c4762a1bSJed Brown    Input Parameters:
532c4762a1bSJed Brown    tao - the Tao context
533c4762a1bSJed Brown    X   - the input vector
534a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
535c4762a1bSJed Brown 
536c4762a1bSJed Brown    Output Parameters:
537c4762a1bSJed Brown    f   - the newly evaluated function
538c4762a1bSJed Brown    G   - the newly evaluated gradient
539c4762a1bSJed Brown */
540c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
541c4762a1bSJed Brown {
542c4762a1bSJed Brown   User              user_ptr = (User)ctx;
543c4762a1bSJed Brown   TS                ts = user_ptr->ts;
544c4762a1bSJed Brown   PetscScalar       *x_ptr,*g;
545c4762a1bSJed Brown   const PetscScalar *y_ptr;
546c4762a1bSJed Brown 
547c4762a1bSJed Brown   PetscFunctionBeginUser;
5485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,&y_ptr));
549c4762a1bSJed Brown   user_ptr->mu = y_ptr[0];
5505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,&y_ptr));
551c4762a1bSJed Brown 
5525f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0.0));
5535f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ts,0));
5545f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,PetscRealConstant(0.001))); /* can be overwritten by command line options */
5555f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
5565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user_ptr->U,&x_ptr));
557c4762a1bSJed Brown   x_ptr[0] = 2.0;
558c4762a1bSJed 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);
5595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user_ptr->U,&x_ptr));
560c4762a1bSJed Brown 
5615f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,user_ptr->U));
562c4762a1bSJed Brown 
5635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(user_ptr->U,&y_ptr));
564c4762a1bSJed 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]);
565c4762a1bSJed Brown 
566c4762a1bSJed Brown   /*   Reset initial conditions for the adjoint integration */
5675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user_ptr->Lambda[0],&x_ptr));
568c4762a1bSJed Brown   x_ptr[0] = 2.*(y_ptr[0]-user_ptr->ob[0]);
569c4762a1bSJed Brown   x_ptr[1] = 2.*(y_ptr[1]-user_ptr->ob[1]);
5705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(user_ptr->U,&y_ptr));
5715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user_ptr->Lambda[0],&x_ptr));
572c4762a1bSJed Brown 
5735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user_ptr->Mup[0],&x_ptr));
574c4762a1bSJed Brown   x_ptr[0] = 0.0;
5755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user_ptr->Mup[0],&x_ptr));
5765f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,1,user_ptr->Lambda,user_ptr->Mup));
577c4762a1bSJed Brown 
5785f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSolve(ts));
579c4762a1bSJed Brown 
5805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user_ptr->Mup[0],&x_ptr));
5815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(user_ptr->Lambda[0],&y_ptr));
5825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(G,&g));
583c4762a1bSJed 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];
5845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user_ptr->Mup[0],&x_ptr));
5855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(user_ptr->Lambda[0],&y_ptr));
5865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(G,&g));
587c4762a1bSJed Brown   PetscFunctionReturn(0);
588c4762a1bSJed Brown }
589c4762a1bSJed Brown 
590c4762a1bSJed Brown PetscErrorCode FormHessian(Tao tao,Vec P,Mat H,Mat Hpre,void *ctx)
591c4762a1bSJed Brown {
592c4762a1bSJed Brown   User           user_ptr = (User)ctx;
593c4762a1bSJed Brown   PetscScalar    harr[1];
594c4762a1bSJed Brown   const PetscInt rows[1] = {0};
595c4762a1bSJed Brown   PetscInt       col = 0;
596c4762a1bSJed Brown 
597c4762a1bSJed Brown   PetscFunctionBeginUser;
5985f80ce2aSJacob Faibussowitsch   CHKERRQ(Adjoint2(P,harr,user_ptr));
5995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(H,1,rows,1,&col,harr,INSERT_VALUES));
600c4762a1bSJed Brown 
6015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
6025f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
603c4762a1bSJed Brown   if (H != Hpre) {
6045f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(Hpre,MAT_FINAL_ASSEMBLY));
6055f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(Hpre,MAT_FINAL_ASSEMBLY));
606c4762a1bSJed Brown   }
607c4762a1bSJed Brown   PetscFunctionReturn(0);
608c4762a1bSJed Brown }
609c4762a1bSJed Brown 
610c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec P,PetscScalar arr[],User ctx)
611c4762a1bSJed Brown {
612c4762a1bSJed Brown   TS                ts = ctx->ts;
613c4762a1bSJed Brown   const PetscScalar *z_ptr;
614c4762a1bSJed Brown   PetscScalar       *x_ptr,*y_ptr,dzdp,dzdp2;
615c4762a1bSJed Brown   Mat               tlmsen;
616c4762a1bSJed Brown 
617c4762a1bSJed Brown   PetscFunctionBeginUser;
618c4762a1bSJed Brown   /* Reset TSAdjoint so that AdjointSetUp will be called again */
6195f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointReset(ts));
620c4762a1bSJed Brown 
621c4762a1bSJed Brown   /* The directional vector should be 1 since it is one-dimensional */
6225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Dir,&x_ptr));
623c4762a1bSJed Brown   x_ptr[0] = 1.;
6245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Dir,&x_ptr));
625c4762a1bSJed Brown 
6265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,&z_ptr));
627c4762a1bSJed Brown   ctx->mu = z_ptr[0];
6285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,&z_ptr));
629c4762a1bSJed Brown 
630c4762a1bSJed Brown   dzdp  = -10.0/(81.0*ctx->mu*ctx->mu) + 2.0*292.0/(2187.0*ctx->mu*ctx->mu*ctx->mu);
631c4762a1bSJed 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);
632c4762a1bSJed Brown 
6335f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0.0));
6345f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ts,0));
6355f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,PetscRealConstant(0.001))); /* can be overwritten by command line options */
6365f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
6375f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostHessianProducts(ts,1,ctx->Lambda2,ctx->Mup2,ctx->Dir));
638c4762a1bSJed Brown 
6395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(ctx->Jacp));
6405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValue(ctx->Jacp,1,0,dzdp,INSERT_VALUES));
6415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(ctx->Jacp,MAT_FINAL_ASSEMBLY));
6425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(ctx->Jacp,MAT_FINAL_ASSEMBLY));
643c4762a1bSJed Brown 
6445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSetForward(ts,ctx->Jacp));
6455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->U,&y_ptr));
646c4762a1bSJed Brown   y_ptr[0] = 2.0;
647c4762a1bSJed Brown   y_ptr[1] = -2.0/3.0 + 10.0/(81.0*ctx->mu) - 292.0/(2187.0*ctx->mu*ctx->mu);
6485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->U,&y_ptr));
6495f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,ctx->U));
650c4762a1bSJed Brown 
651c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
6525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(ctx->U,&z_ptr));
6535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Lambda[0],&y_ptr));
654c4762a1bSJed Brown   y_ptr[0] = 2.*(z_ptr[0]-ctx->ob[0]);
655c4762a1bSJed Brown   y_ptr[1] = 2.*(z_ptr[1]-ctx->ob[1]);
6565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Lambda[0],&y_ptr));
6575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(ctx->U,&z_ptr));
6585f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Mup[0],&y_ptr));
659c4762a1bSJed Brown   y_ptr[0] = 0.0;
6605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Mup[0],&y_ptr));
6615f80ce2aSJacob Faibussowitsch   CHKERRQ(TSForwardGetSensitivities(ts,NULL,&tlmsen));
6625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseGetColumn(tlmsen,0,&x_ptr));
6635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Lambda2[0],&y_ptr));
664c4762a1bSJed Brown   y_ptr[0] = 2.*x_ptr[0];
665c4762a1bSJed Brown   y_ptr[1] = 2.*x_ptr[1];
6665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Lambda2[0],&y_ptr));
6675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Mup2[0],&y_ptr));
668c4762a1bSJed Brown   y_ptr[0] = 0.0;
6695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Mup2[0],&y_ptr));
6705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDenseRestoreColumn(tlmsen,&x_ptr));
6715f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,1,ctx->Lambda,ctx->Mup));
672c4762a1bSJed Brown   if (ctx->implicitform) {
6735f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetIHessianProduct(ts,ctx->Ihp1,IHessianProductUU,ctx->Ihp2,IHessianProductUP,ctx->Ihp3,IHessianProductPU,ctx->Ihp4,IHessianProductPP,ctx));
674c4762a1bSJed Brown   } else {
6755f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSHessianProduct(ts,ctx->Ihp1,RHSHessianProductUU,ctx->Ihp2,RHSHessianProductUP,ctx->Ihp3,RHSHessianProductPU,ctx->Ihp4,RHSHessianProductPP,ctx));
676c4762a1bSJed Brown   }
6775f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSolve(ts));
678c4762a1bSJed Brown 
6795f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Lambda[0],&x_ptr));
6805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->Lambda2[0],&y_ptr));
6815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(ctx->Mup2[0],&z_ptr));
682c4762a1bSJed Brown 
683c4762a1bSJed Brown   arr[0] = x_ptr[1]*dzdp2 + y_ptr[1]*dzdp2 + z_ptr[0];
684c4762a1bSJed Brown 
6855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Lambda2[0],&x_ptr));
6865f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->Lambda2[0],&y_ptr));
6875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(ctx->Mup2[0],&z_ptr));
688c4762a1bSJed Brown 
689c4762a1bSJed Brown   /* Disable second-order adjoint mode */
6905f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointReset(ts));
6915f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointResetForward(ts));
692c4762a1bSJed Brown   PetscFunctionReturn(0);
693c4762a1bSJed Brown }
694c4762a1bSJed Brown 
695c4762a1bSJed Brown /*TEST
696c4762a1bSJed Brown     build:
697c4762a1bSJed Brown       requires: !complex !single
698c4762a1bSJed Brown     test:
699c4762a1bSJed 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
700c4762a1bSJed Brown       output_file: output/ex20opt_p_1.out
701c4762a1bSJed Brown 
702c4762a1bSJed Brown     test:
703c4762a1bSJed Brown       suffix: 2
704c4762a1bSJed 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
705c4762a1bSJed Brown       output_file: output/ex20opt_p_2.out
706c4762a1bSJed Brown 
707c4762a1bSJed Brown     test:
708c4762a1bSJed Brown       suffix: 3
709c4762a1bSJed Brown       args:  -ts_type cn -ts_adapt_type none -mu 100 -ts_dt 0.01 -viewer_binary_skip_info -tao_monitor -tao_view
710c4762a1bSJed Brown       output_file: output/ex20opt_p_3.out
711c4762a1bSJed Brown 
712c4762a1bSJed Brown     test:
713c4762a1bSJed Brown       suffix: 4
714c4762a1bSJed 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
715c4762a1bSJed Brown       output_file: output/ex20opt_p_4.out
716c4762a1bSJed Brown 
717c4762a1bSJed Brown TEST*/
718