xref: /petsc/src/ts/tutorials/optimal_control/ex1.c (revision a82e8c82ed9474375a7f877f23dfa96948657643)
1c4762a1bSJed Brown #include <petsctao.h>
2c4762a1bSJed Brown #include <petscts.h>
3c4762a1bSJed Brown 
4c4762a1bSJed Brown typedef struct _n_aircraft *Aircraft;
5c4762a1bSJed Brown struct _n_aircraft {
6c4762a1bSJed Brown   TS        ts,quadts;
7c4762a1bSJed Brown   Vec       V,W;    /* control variables V and W */
8c4762a1bSJed Brown   PetscInt  nsteps; /* number of time steps */
9c4762a1bSJed Brown   PetscReal ftime;
10c4762a1bSJed Brown   Mat       A,H;
11c4762a1bSJed Brown   Mat       Jacp,DRDU,DRDP;
12c4762a1bSJed Brown   Vec       U,Lambda[1],Mup[1],Lambda2[1],Mup2[1],Dir;
13c4762a1bSJed Brown   Vec       rhshp1[1],rhshp2[1],rhshp3[1],rhshp4[1],inthp1[1],inthp2[1],inthp3[1],inthp4[1];
14c4762a1bSJed Brown   PetscReal lv,lw;
15c4762a1bSJed Brown   PetscBool mf,eh;
16c4762a1bSJed Brown };
17c4762a1bSJed Brown 
18c4762a1bSJed Brown PetscErrorCode FormObjFunctionGradient(Tao,Vec,PetscReal *,Vec,void *);
19c4762a1bSJed Brown PetscErrorCode FormObjHessian(Tao,Vec,Mat,Mat,void *);
20c4762a1bSJed Brown PetscErrorCode ComputeObjHessianWithSOA(Vec,PetscScalar[],Aircraft);
21c4762a1bSJed Brown PetscErrorCode MatrixFreeObjHessian(Tao,Vec,Mat,Mat,void *);
22c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat,Vec,Vec);
23c4762a1bSJed Brown 
24c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
25c4762a1bSJed Brown {
26c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
27c4762a1bSJed Brown   const PetscScalar *u,*v,*w;
28c4762a1bSJed Brown   PetscScalar       *f;
29c4762a1bSJed Brown   PetscInt          step;
30c4762a1bSJed Brown   PetscErrorCode    ierr;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   PetscFunctionBeginUser;
33c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&step);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
35c4762a1bSJed Brown   ierr = VecGetArrayRead(actx->V,&v);CHKERRQ(ierr);
36c4762a1bSJed Brown   ierr = VecGetArrayRead(actx->W,&w);CHKERRQ(ierr);
37c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
38c4762a1bSJed Brown   f[0] = v[step]*PetscCosReal(w[step]);
39c4762a1bSJed Brown   f[1] = v[step]*PetscSinReal(w[step]);
40c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = VecRestoreArrayRead(actx->V,&v);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = VecRestoreArrayRead(actx->W,&w);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
44c4762a1bSJed Brown   PetscFunctionReturn(0);
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
47c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec U,Mat A,void *ctx)
48c4762a1bSJed Brown {
49c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
50c4762a1bSJed Brown   const PetscScalar *u,*v,*w;
51c4762a1bSJed Brown   PetscInt          step,rows[2] = {0,1},rowcol[2];
52c4762a1bSJed Brown   PetscScalar       Jp[2][2];
53c4762a1bSJed Brown   PetscErrorCode    ierr;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   PetscFunctionBeginUser;
56c4762a1bSJed Brown   ierr = MatZeroEntries(A);CHKERRQ(ierr);
57c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&step);CHKERRQ(ierr);
58c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
59c4762a1bSJed Brown   ierr = VecGetArrayRead(actx->V,&v);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = VecGetArrayRead(actx->W,&w);CHKERRQ(ierr);
61c4762a1bSJed Brown 
627b0e2f17SHong Zhang   Jp[0][0] = PetscCosReal(w[step]);
637b0e2f17SHong Zhang   Jp[0][1] = -v[step]*PetscSinReal(w[step]);
647b0e2f17SHong Zhang   Jp[1][0] = PetscSinReal(w[step]);
657b0e2f17SHong Zhang   Jp[1][1] = v[step]*PetscCosReal(w[step]);
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = VecRestoreArrayRead(actx->V,&v);CHKERRQ(ierr);
69c4762a1bSJed Brown   ierr = VecRestoreArrayRead(actx->W,&w);CHKERRQ(ierr);
70c4762a1bSJed Brown 
717b0e2f17SHong Zhang   rowcol[0] = 2*step;
727b0e2f17SHong Zhang   rowcol[1] = 2*step+1;
73c4762a1bSJed Brown   ierr      = MatSetValues(A,2,rows,2,rowcol,&Jp[0][0],INSERT_VALUES);CHKERRQ(ierr);
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
77c4762a1bSJed Brown   PetscFunctionReturn(0);
78c4762a1bSJed Brown }
79c4762a1bSJed Brown 
80c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
81c4762a1bSJed Brown {
82c4762a1bSJed Brown   PetscFunctionBeginUser;
83c4762a1bSJed Brown   PetscFunctionReturn(0);
84c4762a1bSJed Brown }
85c4762a1bSJed Brown 
86c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
87c4762a1bSJed Brown {
88c4762a1bSJed Brown   PetscFunctionBeginUser;
89c4762a1bSJed Brown   PetscFunctionReturn(0);
90c4762a1bSJed Brown }
91c4762a1bSJed Brown 
92c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
93c4762a1bSJed Brown {
94c4762a1bSJed Brown   PetscFunctionBeginUser;
95c4762a1bSJed Brown   PetscFunctionReturn(0);
96c4762a1bSJed Brown }
97c4762a1bSJed Brown 
98c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
99c4762a1bSJed Brown {
100c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
101c4762a1bSJed Brown   const PetscScalar *v,*w,*vl,*vr,*u;
102c4762a1bSJed Brown   PetscScalar       *vhv;
103c4762a1bSJed Brown   PetscScalar       dJpdP[2][2][2]={{{0}}};
104c4762a1bSJed Brown   PetscInt          step,i,j,k;
105c4762a1bSJed Brown   PetscErrorCode    ierr;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   PetscFunctionBeginUser;
108c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&step);CHKERRQ(ierr);
109c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
110c4762a1bSJed Brown   ierr = VecGetArrayRead(actx->V,&v);CHKERRQ(ierr);
111c4762a1bSJed Brown   ierr = VecGetArrayRead(actx->W,&w);CHKERRQ(ierr);
112c4762a1bSJed Brown   ierr = VecGetArrayRead(Vl[0],&vl);CHKERRQ(ierr);
113c4762a1bSJed Brown   ierr = VecGetArrayRead(Vr,&vr);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = VecSet(VHV[0],0.0);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = VecGetArray(VHV[0],&vhv);CHKERRQ(ierr);
116c4762a1bSJed Brown 
1177b0e2f17SHong Zhang   dJpdP[0][0][1] = -PetscSinReal(w[step]);
1187b0e2f17SHong Zhang   dJpdP[0][1][0] = -PetscSinReal(w[step]);
1197b0e2f17SHong Zhang   dJpdP[0][1][1] = -v[step]*PetscCosReal(w[step]);
1207b0e2f17SHong Zhang   dJpdP[1][0][1] = PetscCosReal(w[step]);
1217b0e2f17SHong Zhang   dJpdP[1][1][0] = PetscCosReal(w[step]);
1227b0e2f17SHong Zhang   dJpdP[1][1][1] = -v[step]*PetscSinReal(w[step]);
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   for (j=0; j<2; j++) {
1257b0e2f17SHong Zhang     vhv[2*step+j] = 0;
126c4762a1bSJed Brown     for (k=0; k<2; k++)
127c4762a1bSJed Brown       for (i=0; i<2; i++)
1287b0e2f17SHong Zhang         vhv[2*step+j] += vl[i]*dJpdP[i][j][k]*vr[2*step+k];
129c4762a1bSJed Brown   }
130c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
131c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Vl[0],&vl);CHKERRQ(ierr);
132c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Vr,&vr);CHKERRQ(ierr);
133c4762a1bSJed Brown   ierr = VecRestoreArray(VHV[0],&vhv);CHKERRQ(ierr);
134c4762a1bSJed Brown   PetscFunctionReturn(0);
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown /* Vl in NULL,updates to VHV must be added */
138c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
139c4762a1bSJed Brown {
140c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
141c4762a1bSJed Brown   const PetscScalar *v,*w,*vr,*u;
142c4762a1bSJed Brown   PetscScalar       *vhv;
143c4762a1bSJed Brown   PetscScalar       dRudU[2][2]={{0}};
144c4762a1bSJed Brown   PetscInt          step,j,k;
145c4762a1bSJed Brown   PetscErrorCode    ierr;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   PetscFunctionBeginUser;
148c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&step);CHKERRQ(ierr);
149c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
150c4762a1bSJed Brown   ierr = VecGetArrayRead(actx->V,&v);CHKERRQ(ierr);
151c4762a1bSJed Brown   ierr = VecGetArrayRead(actx->W,&w);CHKERRQ(ierr);
152c4762a1bSJed Brown   ierr = VecGetArrayRead(Vr,&vr);CHKERRQ(ierr);
153c4762a1bSJed Brown   ierr = VecGetArray(VHV[0],&vhv);CHKERRQ(ierr);
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   dRudU[0][0] = 2.0;
156c4762a1bSJed Brown   dRudU[1][1] = 2.0;
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   for (j=0; j<2; j++) {
159c4762a1bSJed Brown     vhv[j] = 0;
160c4762a1bSJed Brown     for (k=0; k<2; k++)
161c4762a1bSJed Brown         vhv[j] += dRudU[j][k]*vr[k];
162c4762a1bSJed Brown   }
163c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
164c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Vr,&vr);CHKERRQ(ierr);
165c4762a1bSJed Brown   ierr = VecRestoreArray(VHV[0],&vhv);CHKERRQ(ierr);
166c4762a1bSJed Brown   PetscFunctionReturn(0);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
170c4762a1bSJed Brown {
171c4762a1bSJed Brown   PetscFunctionBeginUser;
172c4762a1bSJed Brown   PetscFunctionReturn(0);
173c4762a1bSJed Brown }
174c4762a1bSJed Brown 
175c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
176c4762a1bSJed Brown {
177c4762a1bSJed Brown   PetscFunctionBeginUser;
178c4762a1bSJed Brown   PetscFunctionReturn(0);
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
182c4762a1bSJed Brown {
183c4762a1bSJed Brown   PetscFunctionBeginUser;
184c4762a1bSJed Brown   PetscFunctionReturn(0);
185c4762a1bSJed Brown }
186c4762a1bSJed Brown 
187c4762a1bSJed Brown static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,void *ctx)
188c4762a1bSJed Brown {
189c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
190c4762a1bSJed Brown   PetscScalar       *r;
191c4762a1bSJed Brown   PetscReal         dx,dy;
192c4762a1bSJed Brown   const PetscScalar *u;
193c4762a1bSJed Brown   PetscErrorCode    ierr;
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   PetscFunctionBegin;
196c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
197c4762a1bSJed Brown   ierr = VecGetArray(R,&r);CHKERRQ(ierr);
198c4762a1bSJed Brown   dx   = u[0] - actx->lv*t*PetscCosReal(actx->lw);
199c4762a1bSJed Brown   dy   = u[1] - actx->lv*t*PetscSinReal(actx->lw);
200c4762a1bSJed Brown   r[0] = dx*dx+dy*dy;
201c4762a1bSJed Brown   ierr = VecRestoreArray(R,&r);CHKERRQ(ierr);
202c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
203c4762a1bSJed Brown   PetscFunctionReturn(0);
204c4762a1bSJed Brown }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,void *ctx)
207c4762a1bSJed Brown {
208c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
209c4762a1bSJed Brown   PetscScalar       drdu[2][1];
210c4762a1bSJed Brown   const PetscScalar *u;
211c4762a1bSJed Brown   PetscReal         dx,dy;
212c4762a1bSJed Brown   PetscInt          row[] = {0,1},col[] = {0};
213c4762a1bSJed Brown   PetscErrorCode    ierr;
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   PetscFunctionBegin;
216c4762a1bSJed Brown   ierr    = VecGetArrayRead(U,&u);CHKERRQ(ierr);
217c4762a1bSJed Brown   dx      = u[0] - actx->lv*t*PetscCosReal(actx->lw);
218c4762a1bSJed Brown   dy      = u[1] - actx->lv*t*PetscSinReal(actx->lw);
219c4762a1bSJed Brown   drdu[0][0] = 2.*dx;
220c4762a1bSJed Brown   drdu[1][0] = 2.*dy;
221c4762a1bSJed Brown   ierr    = MatSetValues(DRDU,2,row,1,col,&drdu[0][0],INSERT_VALUES);CHKERRQ(ierr);
222c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
223c4762a1bSJed Brown   ierr    = MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
224c4762a1bSJed Brown   ierr    = MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
225c4762a1bSJed Brown   PetscFunctionReturn(0);
226c4762a1bSJed Brown }
227c4762a1bSJed Brown 
228c4762a1bSJed Brown static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,void *ctx)
229c4762a1bSJed Brown {
230c4762a1bSJed Brown   PetscErrorCode ierr;
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   PetscFunctionBegin;
233c4762a1bSJed Brown   ierr = MatZeroEntries(DRDP);CHKERRQ(ierr);
234c4762a1bSJed Brown   ierr = MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
236c4762a1bSJed Brown   PetscFunctionReturn(0);
237c4762a1bSJed Brown }
238c4762a1bSJed Brown 
239c4762a1bSJed Brown int main(int argc,char **argv)
240c4762a1bSJed Brown {
241c4762a1bSJed Brown   Vec                P,PL,PU;
242c4762a1bSJed Brown   struct _n_aircraft aircraft;
243c4762a1bSJed Brown   PetscMPIInt        size;
244c4762a1bSJed Brown   Tao                tao;
245c4762a1bSJed Brown   KSP                ksp;
246c4762a1bSJed Brown   PC                 pc;
247c4762a1bSJed Brown   PetscScalar        *u,*p;
248c4762a1bSJed Brown   PetscInt           i;
249c4762a1bSJed Brown   PetscErrorCode     ierr;
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /* Initialize program */
252c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,NULL);if (ierr) return ierr;
253ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
2542c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* Parameter settings */
257c4762a1bSJed Brown   aircraft.ftime = 1.;   /* time interval in hour */
258c4762a1bSJed Brown   aircraft.nsteps = 10; /* number of steps */
259c4762a1bSJed Brown   aircraft.lv = 2.0; /* leader speed in kmph */
260c4762a1bSJed Brown   aircraft.lw = PETSC_PI/4.; /* leader heading angle */
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   ierr = PetscOptionsGetReal(NULL,NULL,"-ftime",&aircraft.ftime,NULL);CHKERRQ(ierr);
263c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nsteps",&aircraft.nsteps,NULL);CHKERRQ(ierr);
264c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-matrixfree",&aircraft.mf);CHKERRQ(ierr);
265c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-exacthessian",&aircraft.eh);CHKERRQ(ierr);
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
268c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
269c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOBQNLS);CHKERRQ(ierr);
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
272c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&aircraft.A);CHKERRQ(ierr);
273c4762a1bSJed Brown   ierr = MatSetSizes(aircraft.A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
274c4762a1bSJed Brown   ierr = MatSetFromOptions(aircraft.A);CHKERRQ(ierr);
275c4762a1bSJed Brown   ierr = MatSetUp(aircraft.A);CHKERRQ(ierr);
276c4762a1bSJed Brown   ierr = MatAssemblyBegin(aircraft.A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
277c4762a1bSJed Brown   ierr = MatAssemblyEnd(aircraft.A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
278c4762a1bSJed Brown   ierr = MatShift(aircraft.A,1);CHKERRQ(ierr);
279c4762a1bSJed Brown   ierr = MatShift(aircraft.A,-1);CHKERRQ(ierr);
280c4762a1bSJed Brown 
281c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&aircraft.Jacp);CHKERRQ(ierr);
282c4762a1bSJed Brown   ierr = MatSetSizes(aircraft.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,2*aircraft.nsteps);CHKERRQ(ierr);
283c4762a1bSJed Brown   ierr = MatSetFromOptions(aircraft.Jacp);CHKERRQ(ierr);
284c4762a1bSJed Brown   ierr = MatSetUp(aircraft.Jacp);CHKERRQ(ierr);
285c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2*aircraft.nsteps,1,NULL,&aircraft.DRDP);CHKERRQ(ierr);
286c4762a1bSJed Brown   ierr = MatSetUp(aircraft.DRDP);CHKERRQ(ierr);
287c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&aircraft.DRDU);CHKERRQ(ierr);
288c4762a1bSJed Brown   ierr = MatSetUp(aircraft.DRDU);CHKERRQ(ierr);
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   /* Create timestepping solver context */
291c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&aircraft.ts);CHKERRQ(ierr);
292c4762a1bSJed Brown   ierr = TSSetType(aircraft.ts,TSRK);CHKERRQ(ierr);
293c4762a1bSJed Brown   ierr = TSSetRHSFunction(aircraft.ts,NULL,RHSFunction,&aircraft);CHKERRQ(ierr);
294c4762a1bSJed Brown   ierr = TSSetRHSJacobian(aircraft.ts,aircraft.A,aircraft.A,TSComputeRHSJacobianConstant,&aircraft);CHKERRQ(ierr);
295c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(aircraft.ts,aircraft.Jacp,RHSJacobianP,&aircraft);CHKERRQ(ierr);
296c4762a1bSJed Brown   ierr = TSSetExactFinalTime(aircraft.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
297c4762a1bSJed Brown   ierr = TSSetEquationType(aircraft.ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
298c4762a1bSJed Brown 
299c4762a1bSJed Brown   /* Set initial conditions */
300c4762a1bSJed Brown   ierr = MatCreateVecs(aircraft.A,&aircraft.U,NULL);CHKERRQ(ierr);
301c4762a1bSJed Brown   ierr = TSSetSolution(aircraft.ts,aircraft.U);CHKERRQ(ierr);
302c4762a1bSJed Brown   ierr = VecGetArray(aircraft.U,&u);CHKERRQ(ierr);
303c4762a1bSJed Brown   u[0] = 1.5;
304c4762a1bSJed Brown   u[1] = 0;
305c4762a1bSJed Brown   ierr = VecRestoreArray(aircraft.U,&u);CHKERRQ(ierr);
306c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&aircraft.V);CHKERRQ(ierr);
307c4762a1bSJed Brown   ierr = VecSetSizes(aircraft.V,PETSC_DECIDE,aircraft.nsteps);CHKERRQ(ierr);
308c4762a1bSJed Brown   ierr = VecSetUp(aircraft.V);CHKERRQ(ierr);
309c4762a1bSJed Brown   ierr = VecDuplicate(aircraft.V,&aircraft.W);CHKERRQ(ierr);
310c4762a1bSJed Brown   ierr = VecSet(aircraft.V,1.);CHKERRQ(ierr);
311c4762a1bSJed Brown   ierr = VecSet(aircraft.W,PETSC_PI/4.);CHKERRQ(ierr);
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used */
314c4762a1bSJed Brown   ierr = TSSetSaveTrajectory(aircraft.ts);CHKERRQ(ierr);
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   /* Set sensitivity context */
317c4762a1bSJed Brown   ierr = TSCreateQuadratureTS(aircraft.ts,PETSC_FALSE,&aircraft.quadts);CHKERRQ(ierr);
318c4762a1bSJed Brown   ierr = TSSetRHSFunction(aircraft.quadts,NULL,(TSRHSFunction)CostIntegrand,&aircraft);CHKERRQ(ierr);
319c4762a1bSJed Brown   ierr = TSSetRHSJacobian(aircraft.quadts,aircraft.DRDU,aircraft.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&aircraft);CHKERRQ(ierr);
320c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(aircraft.quadts,aircraft.DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&aircraft);CHKERRQ(ierr);
321c4762a1bSJed Brown   ierr = MatCreateVecs(aircraft.A,&aircraft.Lambda[0],NULL);CHKERRQ(ierr);
322c4762a1bSJed Brown   ierr = MatCreateVecs(aircraft.Jacp,&aircraft.Mup[0],NULL);CHKERRQ(ierr);
323c4762a1bSJed Brown   if (aircraft.eh) {
324c4762a1bSJed Brown     ierr = MatCreateVecs(aircraft.A,&aircraft.rhshp1[0],NULL);CHKERRQ(ierr);
325c4762a1bSJed Brown     ierr = MatCreateVecs(aircraft.A,&aircraft.rhshp2[0],NULL);CHKERRQ(ierr);
326c4762a1bSJed Brown     ierr = MatCreateVecs(aircraft.Jacp,&aircraft.rhshp3[0],NULL);CHKERRQ(ierr);
327c4762a1bSJed Brown     ierr = MatCreateVecs(aircraft.Jacp,&aircraft.rhshp4[0],NULL);CHKERRQ(ierr);
328c4762a1bSJed Brown     ierr = MatCreateVecs(aircraft.DRDU,&aircraft.inthp1[0],NULL);CHKERRQ(ierr);
329c4762a1bSJed Brown     ierr = MatCreateVecs(aircraft.DRDU,&aircraft.inthp2[0],NULL);CHKERRQ(ierr);
330c4762a1bSJed Brown     ierr = MatCreateVecs(aircraft.DRDP,&aircraft.inthp3[0],NULL);CHKERRQ(ierr);
331c4762a1bSJed Brown     ierr = MatCreateVecs(aircraft.DRDP,&aircraft.inthp4[0],NULL);CHKERRQ(ierr);
332c4762a1bSJed Brown     ierr = MatCreateVecs(aircraft.Jacp,&aircraft.Dir,NULL);CHKERRQ(ierr);
333c4762a1bSJed Brown     ierr = TSSetRHSHessianProduct(aircraft.ts,aircraft.rhshp1,RHSHessianProductUU,aircraft.rhshp2,RHSHessianProductUP,aircraft.rhshp3,RHSHessianProductPU,aircraft.rhshp4,RHSHessianProductPP,&aircraft);CHKERRQ(ierr);
334c4762a1bSJed Brown     ierr = TSSetRHSHessianProduct(aircraft.quadts,aircraft.inthp1,IntegrandHessianProductUU,aircraft.inthp2,IntegrandHessianProductUP,aircraft.inthp3,IntegrandHessianProductPU,aircraft.inthp4,IntegrandHessianProductPP,&aircraft);CHKERRQ(ierr);
335c4762a1bSJed Brown     ierr = MatCreateVecs(aircraft.A,&aircraft.Lambda2[0],NULL);CHKERRQ(ierr);
336c4762a1bSJed Brown     ierr = MatCreateVecs(aircraft.Jacp,&aircraft.Mup2[0],NULL);CHKERRQ(ierr);
337c4762a1bSJed Brown   }
338c4762a1bSJed Brown   ierr = TSSetFromOptions(aircraft.ts);CHKERRQ(ierr);
339c4762a1bSJed Brown   ierr = TSSetMaxTime(aircraft.ts,aircraft.ftime);CHKERRQ(ierr);
340c4762a1bSJed Brown   ierr = TSSetTimeStep(aircraft.ts,aircraft.ftime/aircraft.nsteps);CHKERRQ(ierr);
341c4762a1bSJed Brown 
342c4762a1bSJed Brown   /* Set initial solution guess */
343c4762a1bSJed Brown   ierr = MatCreateVecs(aircraft.Jacp,&P,NULL);CHKERRQ(ierr);
344c4762a1bSJed Brown   ierr = VecGetArray(P,&p);CHKERRQ(ierr);
345c4762a1bSJed Brown   for (i=0; i<aircraft.nsteps; i++) {
346c4762a1bSJed Brown     p[2*i] = 2.0;
347c4762a1bSJed Brown     p[2*i+1] = PETSC_PI/2.0;
348c4762a1bSJed Brown   }
349c4762a1bSJed Brown   ierr = VecRestoreArray(P,&p);CHKERRQ(ierr);
350c4762a1bSJed Brown   ierr = VecDuplicate(P,&PU);CHKERRQ(ierr);
351c4762a1bSJed Brown   ierr = VecDuplicate(P,&PL);CHKERRQ(ierr);
352c4762a1bSJed Brown   ierr = VecGetArray(PU,&p);CHKERRQ(ierr);
353c4762a1bSJed Brown   for (i=0; i<aircraft.nsteps; i++) {
354c4762a1bSJed Brown     p[2*i] = 2.0;
355c4762a1bSJed Brown     p[2*i+1] = PETSC_PI;
356c4762a1bSJed Brown   }
357c4762a1bSJed Brown   ierr = VecRestoreArray(PU,&p);CHKERRQ(ierr);
358c4762a1bSJed Brown   ierr = VecGetArray(PL,&p);CHKERRQ(ierr);
359c4762a1bSJed Brown   for (i=0; i<aircraft.nsteps; i++) {
360c4762a1bSJed Brown     p[2*i] = 0.0;
361c4762a1bSJed Brown     p[2*i+1] = -PETSC_PI;
362c4762a1bSJed Brown   }
363c4762a1bSJed Brown   ierr = VecRestoreArray(PL,&p);CHKERRQ(ierr);
364c4762a1bSJed Brown 
365*a82e8c82SStefano Zampini   ierr = TaoSetSolution(tao,P);CHKERRQ(ierr);
366c4762a1bSJed Brown   ierr = TaoSetVariableBounds(tao,PL,PU);CHKERRQ(ierr);
367c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
368*a82e8c82SStefano Zampini   ierr = TaoSetObjectiveAndGradient(tao,NULL,FormObjFunctionGradient,(void *)&aircraft);CHKERRQ(ierr);
369c4762a1bSJed Brown 
370c4762a1bSJed Brown   if (aircraft.eh) {
371c4762a1bSJed Brown     if (aircraft.mf) {
372c4762a1bSJed Brown       ierr = MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2*aircraft.nsteps,2*aircraft.nsteps,(void*)&aircraft,&aircraft.H);CHKERRQ(ierr);
373c4762a1bSJed Brown       ierr = MatShellSetOperation(aircraft.H,MATOP_MULT,(void(*)(void))MyMatMult);CHKERRQ(ierr);
374c4762a1bSJed Brown       ierr = MatSetOption(aircraft.H,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
375*a82e8c82SStefano Zampini       ierr = TaoSetHessian(tao,aircraft.H,aircraft.H,MatrixFreeObjHessian,(void*)&aircraft);CHKERRQ(ierr);
376c4762a1bSJed Brown     } else {
377c4762a1bSJed Brown       ierr = MatCreateDense(MPI_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,2*aircraft.nsteps,2*aircraft.nsteps,NULL,&(aircraft.H));CHKERRQ(ierr);
378c4762a1bSJed Brown       ierr = MatSetOption(aircraft.H,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
379*a82e8c82SStefano Zampini       ierr = TaoSetHessian(tao,aircraft.H,aircraft.H,FormObjHessian,(void *)&aircraft);CHKERRQ(ierr);
380c4762a1bSJed Brown     }
381c4762a1bSJed Brown   }
382c4762a1bSJed Brown 
383c4762a1bSJed Brown   /* Check for any TAO command line options */
384c4762a1bSJed Brown   ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
385c4762a1bSJed Brown   if (ksp) {
386c4762a1bSJed Brown     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
387c4762a1bSJed Brown     ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
388c4762a1bSJed Brown   }
389c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
390c4762a1bSJed Brown 
391c4762a1bSJed Brown   ierr = TaoSolve(tao);CHKERRQ(ierr);
392c4762a1bSJed Brown   ierr = VecView(P,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
393c4762a1bSJed Brown 
394c4762a1bSJed Brown   /* Free TAO data structures */
395c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
396c4762a1bSJed Brown 
397c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
398c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
399c4762a1bSJed Brown      are no longer needed.
400c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
401c4762a1bSJed Brown   ierr = TSDestroy(&aircraft.ts);CHKERRQ(ierr);
402c4762a1bSJed Brown   ierr = MatDestroy(&aircraft.A);CHKERRQ(ierr);
403c4762a1bSJed Brown   ierr = VecDestroy(&aircraft.U);CHKERRQ(ierr);
404c4762a1bSJed Brown   ierr = VecDestroy(&aircraft.V);CHKERRQ(ierr);
405c4762a1bSJed Brown   ierr = VecDestroy(&aircraft.W);CHKERRQ(ierr);
406c4762a1bSJed Brown   ierr = VecDestroy(&P);CHKERRQ(ierr);
407c4762a1bSJed Brown   ierr = VecDestroy(&PU);CHKERRQ(ierr);
408c4762a1bSJed Brown   ierr = VecDestroy(&PL);CHKERRQ(ierr);
409c4762a1bSJed Brown   ierr = MatDestroy(&aircraft.Jacp);CHKERRQ(ierr);
410c4762a1bSJed Brown   ierr = MatDestroy(&aircraft.DRDU);CHKERRQ(ierr);
411c4762a1bSJed Brown   ierr = MatDestroy(&aircraft.DRDP);CHKERRQ(ierr);
412c4762a1bSJed Brown   ierr = VecDestroy(&aircraft.Lambda[0]);CHKERRQ(ierr);
413c4762a1bSJed Brown   ierr = VecDestroy(&aircraft.Mup[0]);CHKERRQ(ierr);
414c4762a1bSJed Brown   ierr = VecDestroy(&P);CHKERRQ(ierr);
415c4762a1bSJed Brown   if (aircraft.eh) {
416c4762a1bSJed Brown     ierr = VecDestroy(&aircraft.Lambda2[0]);CHKERRQ(ierr);
417c4762a1bSJed Brown     ierr = VecDestroy(&aircraft.Mup2[0]);CHKERRQ(ierr);
418c4762a1bSJed Brown     ierr = VecDestroy(&aircraft.Dir);CHKERRQ(ierr);
419c4762a1bSJed Brown     ierr = VecDestroy(&aircraft.rhshp1[0]);CHKERRQ(ierr);
420c4762a1bSJed Brown     ierr = VecDestroy(&aircraft.rhshp2[0]);CHKERRQ(ierr);
421c4762a1bSJed Brown     ierr = VecDestroy(&aircraft.rhshp3[0]);CHKERRQ(ierr);
422c4762a1bSJed Brown     ierr = VecDestroy(&aircraft.rhshp4[0]);CHKERRQ(ierr);
423c4762a1bSJed Brown     ierr = VecDestroy(&aircraft.inthp1[0]);CHKERRQ(ierr);
424c4762a1bSJed Brown     ierr = VecDestroy(&aircraft.inthp2[0]);CHKERRQ(ierr);
425c4762a1bSJed Brown     ierr = VecDestroy(&aircraft.inthp3[0]);CHKERRQ(ierr);
426c4762a1bSJed Brown     ierr = VecDestroy(&aircraft.inthp4[0]);CHKERRQ(ierr);
427c4762a1bSJed Brown     ierr = MatDestroy(&aircraft.H);CHKERRQ(ierr);
428c4762a1bSJed Brown   }
429c4762a1bSJed Brown   ierr = PetscFinalize();
430c4762a1bSJed Brown   return ierr;
431c4762a1bSJed Brown }
432c4762a1bSJed Brown 
433c4762a1bSJed Brown /*
434c4762a1bSJed Brown    FormObjFunctionGradient - Evaluates the function and corresponding gradient.
435c4762a1bSJed Brown 
436c4762a1bSJed Brown    Input Parameters:
437c4762a1bSJed Brown    tao - the Tao context
438c4762a1bSJed Brown    P   - the input vector
439*a82e8c82SStefano Zampini    ctx - optional aircraft-defined context, as set by TaoSetObjectiveAndGradient()
440c4762a1bSJed Brown 
441c4762a1bSJed Brown    Output Parameters:
442c4762a1bSJed Brown    f   - the newly evaluated function
443c4762a1bSJed Brown    G   - the newly evaluated gradient
444c4762a1bSJed Brown */
445c4762a1bSJed Brown PetscErrorCode FormObjFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
446c4762a1bSJed Brown {
447c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
448c4762a1bSJed Brown   TS                ts = actx->ts;
449c4762a1bSJed Brown   Vec               Q;
450c4762a1bSJed Brown   const PetscScalar *p,*q;
451c4762a1bSJed Brown   PetscScalar       *u,*v,*w;
452c4762a1bSJed Brown   PetscInt          i;
453c4762a1bSJed Brown   PetscErrorCode    ierr;
454c4762a1bSJed Brown 
455c4762a1bSJed Brown   PetscFunctionBeginUser;
456c4762a1bSJed Brown   ierr = VecGetArrayRead(P,&p);CHKERRQ(ierr);
457c4762a1bSJed Brown   ierr = VecGetArray(actx->V,&v);CHKERRQ(ierr);
458c4762a1bSJed Brown   ierr = VecGetArray(actx->W,&w);CHKERRQ(ierr);
459c4762a1bSJed Brown   for (i=0; i<actx->nsteps; i++) {
460c4762a1bSJed Brown     v[i] = p[2*i];
461c4762a1bSJed Brown     w[i] = p[2*i+1];
462c4762a1bSJed Brown   }
463c4762a1bSJed Brown   ierr = VecRestoreArrayRead(P,&p);CHKERRQ(ierr);
464c4762a1bSJed Brown   ierr = VecRestoreArray(actx->V,&v);CHKERRQ(ierr);
465c4762a1bSJed Brown   ierr = VecRestoreArray(actx->W,&w);CHKERRQ(ierr);
466c4762a1bSJed Brown 
467c4762a1bSJed Brown   ierr = TSSetTime(ts,0.0);CHKERRQ(ierr);
468c4762a1bSJed Brown   ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr);
469c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
470c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,actx->ftime/actx->nsteps);CHKERRQ(ierr);
471c4762a1bSJed Brown 
472c4762a1bSJed Brown   /* reinitialize system state */
473c4762a1bSJed Brown   ierr = VecGetArray(actx->U,&u);CHKERRQ(ierr);
474c4762a1bSJed Brown   u[0] = 2.0;
475c4762a1bSJed Brown   u[1] = 0;
476c4762a1bSJed Brown   ierr = VecRestoreArray(actx->U,&u);CHKERRQ(ierr);
477c4762a1bSJed Brown 
478c4762a1bSJed Brown   /* reinitialize the integral value */
479c4762a1bSJed Brown   ierr = TSGetCostIntegral(ts,&Q);CHKERRQ(ierr);
480c4762a1bSJed Brown   ierr = VecSet(Q,0.0);CHKERRQ(ierr);
481c4762a1bSJed Brown 
482c4762a1bSJed Brown   ierr = TSSolve(ts,actx->U);CHKERRQ(ierr);
483c4762a1bSJed Brown 
484c4762a1bSJed Brown   /* Reset initial conditions for the adjoint integration */
485c4762a1bSJed Brown   ierr = VecSet(actx->Lambda[0],0.0);CHKERRQ(ierr);
486c4762a1bSJed Brown   ierr = VecSet(actx->Mup[0],0.0);CHKERRQ(ierr);
487c4762a1bSJed Brown   ierr = TSSetCostGradients(ts,1,actx->Lambda,actx->Mup);CHKERRQ(ierr);
488c4762a1bSJed Brown 
489c4762a1bSJed Brown   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
490c4762a1bSJed Brown   ierr = VecCopy(actx->Mup[0],G);CHKERRQ(ierr);
491c4762a1bSJed Brown   ierr = TSGetCostIntegral(ts,&Q);CHKERRQ(ierr);
492c4762a1bSJed Brown   ierr = VecGetArrayRead(Q,&q);CHKERRQ(ierr);
493c4762a1bSJed Brown   *f   = q[0];
494c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Q,&q);CHKERRQ(ierr);
495c4762a1bSJed Brown   PetscFunctionReturn(0);
496c4762a1bSJed Brown }
497c4762a1bSJed Brown 
498c4762a1bSJed Brown PetscErrorCode FormObjHessian(Tao tao,Vec P,Mat H,Mat Hpre,void *ctx)
499c4762a1bSJed Brown {
500c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
501c4762a1bSJed Brown   const PetscScalar *p;
502c4762a1bSJed Brown   PetscScalar       *harr,*v,*w,one = 1.0;
503c4762a1bSJed Brown   PetscInt          ind[1];
504c4762a1bSJed Brown   PetscInt          *cols,i;
505c4762a1bSJed Brown   Vec               Dir;
506c4762a1bSJed Brown   PetscErrorCode    ierr;
507c4762a1bSJed Brown 
508c4762a1bSJed Brown   PetscFunctionBeginUser;
509c4762a1bSJed Brown   /* set up control parameters */
510c4762a1bSJed Brown   ierr = VecGetArrayRead(P,&p);CHKERRQ(ierr);
511c4762a1bSJed Brown   ierr = VecGetArray(actx->V,&v);CHKERRQ(ierr);
512c4762a1bSJed Brown   ierr = VecGetArray(actx->W,&w);CHKERRQ(ierr);
513c4762a1bSJed Brown   for (i=0; i<actx->nsteps; i++) {
514c4762a1bSJed Brown     v[i] = p[2*i];
515c4762a1bSJed Brown     w[i] = p[2*i+1];
516c4762a1bSJed Brown   }
517c4762a1bSJed Brown   ierr = VecRestoreArrayRead(P,&p);CHKERRQ(ierr);
518c4762a1bSJed Brown   ierr = VecRestoreArray(actx->V,&v);CHKERRQ(ierr);
519c4762a1bSJed Brown   ierr = VecRestoreArray(actx->W,&w);CHKERRQ(ierr);
520c4762a1bSJed Brown 
521c4762a1bSJed Brown   ierr = PetscMalloc1(2*actx->nsteps,&harr);CHKERRQ(ierr);
522c4762a1bSJed Brown   ierr = PetscMalloc1(2*actx->nsteps,&cols);CHKERRQ(ierr);
523c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) cols[i] = i;
524c4762a1bSJed Brown   ierr = VecDuplicate(P,&Dir);CHKERRQ(ierr);
525c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) {
526c4762a1bSJed Brown     ind[0] = i;
527c4762a1bSJed Brown     ierr   = VecSet(Dir,0.0);CHKERRQ(ierr);
528c4762a1bSJed Brown     ierr   = VecSetValues(Dir,1,ind,&one,INSERT_VALUES);CHKERRQ(ierr);
529c4762a1bSJed Brown     ierr   = VecAssemblyBegin(Dir);CHKERRQ(ierr);
530c4762a1bSJed Brown     ierr   = VecAssemblyEnd(Dir);CHKERRQ(ierr);
531c4762a1bSJed Brown     ierr   = ComputeObjHessianWithSOA(Dir,harr,actx);CHKERRQ(ierr);
532c4762a1bSJed Brown     ierr   = MatSetValues(H,1,ind,2*actx->nsteps,cols,harr,INSERT_VALUES);CHKERRQ(ierr);
533c4762a1bSJed Brown     ierr   = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
534c4762a1bSJed Brown     ierr   = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
535c4762a1bSJed Brown     if (H != Hpre) {
536c4762a1bSJed Brown       ierr   = MatAssemblyBegin(Hpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
537c4762a1bSJed Brown       ierr   = MatAssemblyEnd(Hpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
538c4762a1bSJed Brown     }
539c4762a1bSJed Brown   }
540c4762a1bSJed Brown   ierr = PetscFree(cols);CHKERRQ(ierr);
541c4762a1bSJed Brown   ierr = PetscFree(harr);CHKERRQ(ierr);
542c4762a1bSJed Brown   ierr = VecDestroy(&Dir);CHKERRQ(ierr);
543c4762a1bSJed Brown   PetscFunctionReturn(0);
544c4762a1bSJed Brown }
545c4762a1bSJed Brown 
546c4762a1bSJed Brown PetscErrorCode MatrixFreeObjHessian(Tao tao, Vec P, Mat H, Mat Hpre, void *ctx)
547c4762a1bSJed Brown {
548c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
549c4762a1bSJed Brown   PetscScalar       *v,*w;
550c4762a1bSJed Brown   const PetscScalar *p;
551c4762a1bSJed Brown   PetscInt          i;
552c4762a1bSJed Brown   PetscErrorCode    ierr;
553c4762a1bSJed Brown 
554c4762a1bSJed Brown   PetscFunctionBegin;
555c4762a1bSJed Brown   ierr = VecGetArrayRead(P,&p);CHKERRQ(ierr);
556c4762a1bSJed Brown   ierr = VecGetArray(actx->V,&v);CHKERRQ(ierr);
557c4762a1bSJed Brown   ierr = VecGetArray(actx->W,&w);CHKERRQ(ierr);
558c4762a1bSJed Brown   for (i=0; i<actx->nsteps; i++) {
559c4762a1bSJed Brown     v[i] = p[2*i];
560c4762a1bSJed Brown     w[i] = p[2*i+1];
561c4762a1bSJed Brown   }
562c4762a1bSJed Brown   ierr = VecRestoreArrayRead(P,&p);CHKERRQ(ierr);
563c4762a1bSJed Brown   ierr = VecRestoreArray(actx->V,&v);CHKERRQ(ierr);
564c4762a1bSJed Brown   ierr = VecRestoreArray(actx->W,&w);CHKERRQ(ierr);
565c4762a1bSJed Brown   PetscFunctionReturn(0);
566c4762a1bSJed Brown }
567c4762a1bSJed Brown 
568c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y)
569c4762a1bSJed Brown {
570c4762a1bSJed Brown   PetscScalar    *y;
571c4762a1bSJed Brown   void           *ptr;
572c4762a1bSJed Brown   PetscErrorCode ierr;
573c4762a1bSJed Brown 
574c4762a1bSJed Brown   PetscFunctionBegin;
575c4762a1bSJed Brown   ierr = MatShellGetContext(H_shell,&ptr);CHKERRQ(ierr);
576c4762a1bSJed Brown   ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
577c4762a1bSJed Brown   ierr = ComputeObjHessianWithSOA(X,y,(Aircraft)ptr);CHKERRQ(ierr);
578c4762a1bSJed Brown   ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
579c4762a1bSJed Brown   PetscFunctionReturn(0);
580c4762a1bSJed Brown }
581c4762a1bSJed Brown 
582c4762a1bSJed Brown PetscErrorCode ComputeObjHessianWithSOA(Vec Dir,PetscScalar arr[],Aircraft actx)
583c4762a1bSJed Brown {
584c4762a1bSJed Brown   TS                ts = actx->ts;
585c4762a1bSJed Brown   const PetscScalar *z_ptr;
586c4762a1bSJed Brown   PetscScalar       *u;
587c4762a1bSJed Brown   Vec               Q;
588c4762a1bSJed Brown   PetscInt          i;
589c4762a1bSJed Brown   PetscErrorCode    ierr;
590c4762a1bSJed Brown 
591c4762a1bSJed Brown   PetscFunctionBeginUser;
592c4762a1bSJed Brown   /* Reset TSAdjoint so that AdjointSetUp will be called again */
593c4762a1bSJed Brown   ierr = TSAdjointReset(ts);CHKERRQ(ierr);
594c4762a1bSJed Brown 
595c4762a1bSJed Brown   ierr = TSSetTime(ts,0.0);CHKERRQ(ierr);
596c4762a1bSJed Brown   ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr);
597c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
598c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,actx->ftime/actx->nsteps);CHKERRQ(ierr);
599c4762a1bSJed Brown   ierr = TSSetCostHessianProducts(actx->ts,1,actx->Lambda2,actx->Mup2,Dir);CHKERRQ(ierr);
600c4762a1bSJed Brown 
601c4762a1bSJed Brown   /* reinitialize system state */
602c4762a1bSJed Brown   ierr = VecGetArray(actx->U,&u);CHKERRQ(ierr);
603c4762a1bSJed Brown   u[0] = 2.0;
604c4762a1bSJed Brown   u[1] = 0;
605c4762a1bSJed Brown   ierr = VecRestoreArray(actx->U,&u);CHKERRQ(ierr);
606c4762a1bSJed Brown 
607c4762a1bSJed Brown   /* reinitialize the integral value */
608c4762a1bSJed Brown   ierr = TSGetCostIntegral(ts,&Q);CHKERRQ(ierr);
609c4762a1bSJed Brown   ierr = VecSet(Q,0.0);CHKERRQ(ierr);
610c4762a1bSJed Brown 
611c4762a1bSJed Brown   /* initialize tlm variable */
612c4762a1bSJed Brown   ierr = MatZeroEntries(actx->Jacp);CHKERRQ(ierr);
613c4762a1bSJed Brown   ierr = MatAssemblyBegin(actx->Jacp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
614c4762a1bSJed Brown   ierr = MatAssemblyEnd(actx->Jacp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
615c4762a1bSJed Brown   ierr = TSAdjointSetForward(ts,actx->Jacp);CHKERRQ(ierr);
616c4762a1bSJed Brown 
617c4762a1bSJed Brown   ierr = TSSolve(ts,actx->U);CHKERRQ(ierr);
618c4762a1bSJed Brown 
619c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
620c4762a1bSJed Brown   ierr = VecSet(actx->Lambda[0],0.0);CHKERRQ(ierr);
621c4762a1bSJed Brown   ierr = VecSet(actx->Mup[0],0.0);CHKERRQ(ierr);
622c4762a1bSJed Brown   ierr = VecSet(actx->Lambda2[0],0.0);CHKERRQ(ierr);
623c4762a1bSJed Brown   ierr = VecSet(actx->Mup2[0],0.0);CHKERRQ(ierr);
624c4762a1bSJed Brown   ierr = TSSetCostGradients(ts,1,actx->Lambda,actx->Mup);CHKERRQ(ierr);
625c4762a1bSJed Brown 
626c4762a1bSJed Brown   ierr = TSGetCostIntegral(ts,&Q);CHKERRQ(ierr);
627c4762a1bSJed Brown 
628c4762a1bSJed Brown   /* Reset initial conditions for the adjoint integration */
629c4762a1bSJed Brown   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
630c4762a1bSJed Brown 
631c4762a1bSJed Brown   /* initial condition does not depend on p, so that lambda is not needed to assemble G */
632c4762a1bSJed Brown   ierr = VecGetArrayRead(actx->Mup2[0],&z_ptr);CHKERRQ(ierr);
633c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) arr[i] = z_ptr[i];
634c4762a1bSJed Brown   ierr = VecRestoreArrayRead(actx->Mup2[0],&z_ptr);CHKERRQ(ierr);
635c4762a1bSJed Brown 
636c4762a1bSJed Brown   /* Disable second-order adjoint mode */
637c4762a1bSJed Brown   ierr = TSAdjointReset(ts);CHKERRQ(ierr);
638c4762a1bSJed Brown   ierr = TSAdjointResetForward(ts);CHKERRQ(ierr);
639c4762a1bSJed Brown   PetscFunctionReturn(0);
640c4762a1bSJed Brown }
641c4762a1bSJed Brown 
642c4762a1bSJed Brown /*TEST
643c4762a1bSJed Brown     build:
644c4762a1bSJed Brown       requires: !complex !single
645c4762a1bSJed Brown 
646c4762a1bSJed Brown     test:
647c4762a1bSJed Brown       args:  -ts_adapt_type none -ts_type rk -ts_rk_type 3 -viewer_binary_skip_info -tao_monitor -tao_gatol 1e-7
648c4762a1bSJed Brown 
649c4762a1bSJed Brown     test:
650c4762a1bSJed Brown       suffix: 2
651c4762a1bSJed Brown       args:  -ts_adapt_type none -ts_type rk -ts_rk_type 3 -viewer_binary_skip_info -tao_monitor -tao_view -tao_type bntr -tao_bnk_pc_type none -exacthessian
652c4762a1bSJed Brown 
653c4762a1bSJed Brown     test:
654c4762a1bSJed Brown       suffix: 3
655c4762a1bSJed Brown       args:  -ts_adapt_type none -ts_type rk -ts_rk_type 3 -viewer_binary_skip_info -tao_monitor -tao_view -tao_type bntr -tao_bnk_pc_type none -exacthessian -matrixfree
656c4762a1bSJed Brown TEST*/
657c4762a1bSJed Brown 
658