xref: /petsc/src/ts/tutorials/optimal_control/ex1.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
31c4762a1bSJed Brown   PetscFunctionBeginUser;
329566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&step));
339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->V,&v));
359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->W,&w));
369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
37c4762a1bSJed Brown   f[0] = v[step]*PetscCosReal(w[step]);
38c4762a1bSJed Brown   f[1] = v[step]*PetscSinReal(w[step]);
399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(actx->V,&v));
419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(actx->W,&w));
429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
43c4762a1bSJed Brown   PetscFunctionReturn(0);
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec U,Mat A,void *ctx)
47c4762a1bSJed Brown {
48c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
49c4762a1bSJed Brown   const PetscScalar *u,*v,*w;
50c4762a1bSJed Brown   PetscInt          step,rows[2] = {0,1},rowcol[2];
51c4762a1bSJed Brown   PetscScalar       Jp[2][2];
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   PetscFunctionBeginUser;
549566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(A));
559566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&step));
569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->V,&v));
589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->W,&w));
59c4762a1bSJed Brown 
607b0e2f17SHong Zhang   Jp[0][0] = PetscCosReal(w[step]);
617b0e2f17SHong Zhang   Jp[0][1] = -v[step]*PetscSinReal(w[step]);
627b0e2f17SHong Zhang   Jp[1][0] = PetscSinReal(w[step]);
637b0e2f17SHong Zhang   Jp[1][1] = v[step]*PetscCosReal(w[step]);
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(actx->V,&v));
679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(actx->W,&w));
68c4762a1bSJed Brown 
697b0e2f17SHong Zhang   rowcol[0] = 2*step;
707b0e2f17SHong Zhang   rowcol[1] = 2*step+1;
719566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,2,rows,2,rowcol,&Jp[0][0],INSERT_VALUES));
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
75c4762a1bSJed Brown   PetscFunctionReturn(0);
76c4762a1bSJed Brown }
77c4762a1bSJed Brown 
78c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
79c4762a1bSJed Brown {
80c4762a1bSJed Brown   PetscFunctionBeginUser;
81c4762a1bSJed Brown   PetscFunctionReturn(0);
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
84c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
85c4762a1bSJed Brown {
86c4762a1bSJed Brown   PetscFunctionBeginUser;
87c4762a1bSJed Brown   PetscFunctionReturn(0);
88c4762a1bSJed Brown }
89c4762a1bSJed Brown 
90c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
91c4762a1bSJed Brown {
92c4762a1bSJed Brown   PetscFunctionBeginUser;
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
97c4762a1bSJed Brown {
98c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
99c4762a1bSJed Brown   const PetscScalar *v,*w,*vl,*vr,*u;
100c4762a1bSJed Brown   PetscScalar       *vhv;
101c4762a1bSJed Brown   PetscScalar       dJpdP[2][2][2]={{{0}}};
102c4762a1bSJed Brown   PetscInt          step,i,j,k;
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   PetscFunctionBeginUser;
1059566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&step));
1069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1079566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->V,&v));
1089566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->W,&w));
1099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0],&vl));
1109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr,&vr));
1119566063dSJacob Faibussowitsch   PetscCall(VecSet(VHV[0],0.0));
1129566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0],&vhv));
113c4762a1bSJed Brown 
1147b0e2f17SHong Zhang   dJpdP[0][0][1] = -PetscSinReal(w[step]);
1157b0e2f17SHong Zhang   dJpdP[0][1][0] = -PetscSinReal(w[step]);
1167b0e2f17SHong Zhang   dJpdP[0][1][1] = -v[step]*PetscCosReal(w[step]);
1177b0e2f17SHong Zhang   dJpdP[1][0][1] = PetscCosReal(w[step]);
1187b0e2f17SHong Zhang   dJpdP[1][1][0] = PetscCosReal(w[step]);
1197b0e2f17SHong Zhang   dJpdP[1][1][1] = -v[step]*PetscSinReal(w[step]);
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   for (j=0; j<2; j++) {
1227b0e2f17SHong Zhang     vhv[2*step+j] = 0;
123c4762a1bSJed Brown     for (k=0; k<2; k++)
124c4762a1bSJed Brown       for (i=0; i<2; i++)
1257b0e2f17SHong Zhang         vhv[2*step+j] += vl[i]*dJpdP[i][j][k]*vr[2*step+k];
126c4762a1bSJed Brown   }
1279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0],&vl));
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr,&vr));
1309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0],&vhv));
131c4762a1bSJed Brown   PetscFunctionReturn(0);
132c4762a1bSJed Brown }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown /* Vl in NULL,updates to VHV must be added */
135c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
136c4762a1bSJed Brown {
137c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
138c4762a1bSJed Brown   const PetscScalar *v,*w,*vr,*u;
139c4762a1bSJed Brown   PetscScalar       *vhv;
140c4762a1bSJed Brown   PetscScalar       dRudU[2][2]={{0}};
141c4762a1bSJed Brown   PetscInt          step,j,k;
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   PetscFunctionBeginUser;
1449566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&step));
1459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1469566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->V,&v));
1479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->W,&w));
1489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr,&vr));
1499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0],&vhv));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   dRudU[0][0] = 2.0;
152c4762a1bSJed Brown   dRudU[1][1] = 2.0;
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   for (j=0; j<2; j++) {
155c4762a1bSJed Brown     vhv[j] = 0;
156c4762a1bSJed Brown     for (k=0; k<2; k++)
157c4762a1bSJed Brown         vhv[j] += dRudU[j][k]*vr[k];
158c4762a1bSJed Brown   }
1599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr,&vr));
1619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0],&vhv));
162c4762a1bSJed Brown   PetscFunctionReturn(0);
163c4762a1bSJed Brown }
164c4762a1bSJed Brown 
165c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
166c4762a1bSJed Brown {
167c4762a1bSJed Brown   PetscFunctionBeginUser;
168c4762a1bSJed Brown   PetscFunctionReturn(0);
169c4762a1bSJed Brown }
170c4762a1bSJed Brown 
171c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
172c4762a1bSJed Brown {
173c4762a1bSJed Brown   PetscFunctionBeginUser;
174c4762a1bSJed Brown   PetscFunctionReturn(0);
175c4762a1bSJed Brown }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx)
178c4762a1bSJed Brown {
179c4762a1bSJed Brown   PetscFunctionBeginUser;
180c4762a1bSJed Brown   PetscFunctionReturn(0);
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
183c4762a1bSJed Brown static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,void *ctx)
184c4762a1bSJed Brown {
185c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
186c4762a1bSJed Brown   PetscScalar       *r;
187c4762a1bSJed Brown   PetscReal         dx,dy;
188c4762a1bSJed Brown   const PetscScalar *u;
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   PetscFunctionBegin;
1919566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1929566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R,&r));
193c4762a1bSJed Brown   dx   = u[0] - actx->lv*t*PetscCosReal(actx->lw);
194c4762a1bSJed Brown   dy   = u[1] - actx->lv*t*PetscSinReal(actx->lw);
195c4762a1bSJed Brown   r[0] = dx*dx+dy*dy;
1969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R,&r));
1979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
198c4762a1bSJed Brown   PetscFunctionReturn(0);
199c4762a1bSJed Brown }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,void *ctx)
202c4762a1bSJed Brown {
203c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
204c4762a1bSJed Brown   PetscScalar       drdu[2][1];
205c4762a1bSJed Brown   const PetscScalar *u;
206c4762a1bSJed Brown   PetscReal         dx,dy;
207c4762a1bSJed Brown   PetscInt          row[] = {0,1},col[] = {0};
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   PetscFunctionBegin;
2109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
211c4762a1bSJed Brown   dx      = u[0] - actx->lv*t*PetscCosReal(actx->lw);
212c4762a1bSJed Brown   dy      = u[1] - actx->lv*t*PetscSinReal(actx->lw);
213c4762a1bSJed Brown   drdu[0][0] = 2.*dx;
214c4762a1bSJed Brown   drdu[1][0] = 2.*dy;
2159566063dSJacob Faibussowitsch   PetscCall(MatSetValues(DRDU,2,row,1,col,&drdu[0][0],INSERT_VALUES));
2169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
2179566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY));
2189566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY));
219c4762a1bSJed Brown   PetscFunctionReturn(0);
220c4762a1bSJed Brown }
221c4762a1bSJed Brown 
222c4762a1bSJed Brown static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,void *ctx)
223c4762a1bSJed Brown {
224c4762a1bSJed Brown   PetscFunctionBegin;
2259566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(DRDP));
2269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY));
2279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY));
228c4762a1bSJed Brown   PetscFunctionReturn(0);
229c4762a1bSJed Brown }
230c4762a1bSJed Brown 
231c4762a1bSJed Brown int main(int argc,char **argv)
232c4762a1bSJed Brown {
233c4762a1bSJed Brown   Vec                P,PL,PU;
234c4762a1bSJed Brown   struct _n_aircraft aircraft;
235c4762a1bSJed Brown   PetscMPIInt        size;
236c4762a1bSJed Brown   Tao                tao;
237c4762a1bSJed Brown   KSP                ksp;
238c4762a1bSJed Brown   PC                 pc;
239c4762a1bSJed Brown   PetscScalar        *u,*p;
240c4762a1bSJed Brown   PetscInt           i;
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /* Initialize program */
243*327415f7SBarry Smith   PetscFunctionBeginUser;
2449566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,NULL));
2459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
2463c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
247c4762a1bSJed Brown 
248c4762a1bSJed Brown   /* Parameter settings */
249c4762a1bSJed Brown   aircraft.ftime = 1.;   /* time interval in hour */
250c4762a1bSJed Brown   aircraft.nsteps = 10; /* number of steps */
251c4762a1bSJed Brown   aircraft.lv = 2.0; /* leader speed in kmph */
252c4762a1bSJed Brown   aircraft.lw = PETSC_PI/4.; /* leader heading angle */
253c4762a1bSJed Brown 
2549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-ftime",&aircraft.ftime,NULL));
2559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nsteps",&aircraft.nsteps,NULL));
2569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-matrixfree",&aircraft.mf));
2579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-exacthessian",&aircraft.eh));
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
2609566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
2619566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOBQNLS));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
2649566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&aircraft.A));
2659566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(aircraft.A,PETSC_DECIDE,PETSC_DECIDE,2,2));
2669566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(aircraft.A));
2679566063dSJacob Faibussowitsch   PetscCall(MatSetUp(aircraft.A));
2689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(aircraft.A,MAT_FINAL_ASSEMBLY));
2699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(aircraft.A,MAT_FINAL_ASSEMBLY));
2709566063dSJacob Faibussowitsch   PetscCall(MatShift(aircraft.A,1));
2719566063dSJacob Faibussowitsch   PetscCall(MatShift(aircraft.A,-1));
272c4762a1bSJed Brown 
2739566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&aircraft.Jacp));
2749566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(aircraft.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,2*aircraft.nsteps));
2759566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(aircraft.Jacp));
2769566063dSJacob Faibussowitsch   PetscCall(MatSetUp(aircraft.Jacp));
2779566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2*aircraft.nsteps,1,NULL,&aircraft.DRDP));
2789566063dSJacob Faibussowitsch   PetscCall(MatSetUp(aircraft.DRDP));
2799566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&aircraft.DRDU));
2809566063dSJacob Faibussowitsch   PetscCall(MatSetUp(aircraft.DRDU));
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /* Create timestepping solver context */
2839566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&aircraft.ts));
2849566063dSJacob Faibussowitsch   PetscCall(TSSetType(aircraft.ts,TSRK));
2859566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(aircraft.ts,NULL,RHSFunction,&aircraft));
2869566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(aircraft.ts,aircraft.A,aircraft.A,TSComputeRHSJacobianConstant,&aircraft));
2879566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(aircraft.ts,aircraft.Jacp,RHSJacobianP,&aircraft));
2889566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(aircraft.ts,TS_EXACTFINALTIME_MATCHSTEP));
2899566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(aircraft.ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
290c4762a1bSJed Brown 
291c4762a1bSJed Brown   /* Set initial conditions */
2929566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(aircraft.A,&aircraft.U,NULL));
2939566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(aircraft.ts,aircraft.U));
2949566063dSJacob Faibussowitsch   PetscCall(VecGetArray(aircraft.U,&u));
295c4762a1bSJed Brown   u[0] = 1.5;
296c4762a1bSJed Brown   u[1] = 0;
2979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(aircraft.U,&u));
2989566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&aircraft.V));
2999566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(aircraft.V,PETSC_DECIDE,aircraft.nsteps));
3009566063dSJacob Faibussowitsch   PetscCall(VecSetUp(aircraft.V));
3019566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(aircraft.V,&aircraft.W));
3029566063dSJacob Faibussowitsch   PetscCall(VecSet(aircraft.V,1.));
3039566063dSJacob Faibussowitsch   PetscCall(VecSet(aircraft.W,PETSC_PI/4.));
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used */
3069566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(aircraft.ts));
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   /* Set sensitivity context */
3099566063dSJacob Faibussowitsch   PetscCall(TSCreateQuadratureTS(aircraft.ts,PETSC_FALSE,&aircraft.quadts));
3109566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(aircraft.quadts,NULL,(TSRHSFunction)CostIntegrand,&aircraft));
3119566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(aircraft.quadts,aircraft.DRDU,aircraft.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&aircraft));
3129566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(aircraft.quadts,aircraft.DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&aircraft));
3139566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(aircraft.A,&aircraft.Lambda[0],NULL));
3149566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(aircraft.Jacp,&aircraft.Mup[0],NULL));
315c4762a1bSJed Brown   if (aircraft.eh) {
3169566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.A,&aircraft.rhshp1[0],NULL));
3179566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.A,&aircraft.rhshp2[0],NULL));
3189566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.Jacp,&aircraft.rhshp3[0],NULL));
3199566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.Jacp,&aircraft.rhshp4[0],NULL));
3209566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.DRDU,&aircraft.inthp1[0],NULL));
3219566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.DRDU,&aircraft.inthp2[0],NULL));
3229566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.DRDP,&aircraft.inthp3[0],NULL));
3239566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.DRDP,&aircraft.inthp4[0],NULL));
3249566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.Jacp,&aircraft.Dir,NULL));
3259566063dSJacob Faibussowitsch     PetscCall(TSSetRHSHessianProduct(aircraft.ts,aircraft.rhshp1,RHSHessianProductUU,aircraft.rhshp2,RHSHessianProductUP,aircraft.rhshp3,RHSHessianProductPU,aircraft.rhshp4,RHSHessianProductPP,&aircraft));
3269566063dSJacob Faibussowitsch     PetscCall(TSSetRHSHessianProduct(aircraft.quadts,aircraft.inthp1,IntegrandHessianProductUU,aircraft.inthp2,IntegrandHessianProductUP,aircraft.inthp3,IntegrandHessianProductPU,aircraft.inthp4,IntegrandHessianProductPP,&aircraft));
3279566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.A,&aircraft.Lambda2[0],NULL));
3289566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.Jacp,&aircraft.Mup2[0],NULL));
329c4762a1bSJed Brown   }
3309566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(aircraft.ts));
3319566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(aircraft.ts,aircraft.ftime));
3329566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(aircraft.ts,aircraft.ftime/aircraft.nsteps));
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   /* Set initial solution guess */
3359566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(aircraft.Jacp,&P,NULL));
3369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(P,&p));
337c4762a1bSJed Brown   for (i=0; i<aircraft.nsteps; i++) {
338c4762a1bSJed Brown     p[2*i] = 2.0;
339c4762a1bSJed Brown     p[2*i+1] = PETSC_PI/2.0;
340c4762a1bSJed Brown   }
3419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(P,&p));
3429566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(P,&PU));
3439566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(P,&PL));
3449566063dSJacob Faibussowitsch   PetscCall(VecGetArray(PU,&p));
345c4762a1bSJed Brown   for (i=0; i<aircraft.nsteps; i++) {
346c4762a1bSJed Brown     p[2*i] = 2.0;
347c4762a1bSJed Brown     p[2*i+1] = PETSC_PI;
348c4762a1bSJed Brown   }
3499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(PU,&p));
3509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(PL,&p));
351c4762a1bSJed Brown   for (i=0; i<aircraft.nsteps; i++) {
352c4762a1bSJed Brown     p[2*i] = 0.0;
353c4762a1bSJed Brown     p[2*i+1] = -PETSC_PI;
354c4762a1bSJed Brown   }
3559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(PL,&p));
356c4762a1bSJed Brown 
3579566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,P));
3589566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao,PL,PU));
359c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
3609566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormObjFunctionGradient,(void *)&aircraft));
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   if (aircraft.eh) {
363c4762a1bSJed Brown     if (aircraft.mf) {
3649566063dSJacob Faibussowitsch       PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2*aircraft.nsteps,2*aircraft.nsteps,(void*)&aircraft,&aircraft.H));
3659566063dSJacob Faibussowitsch       PetscCall(MatShellSetOperation(aircraft.H,MATOP_MULT,(void(*)(void))MyMatMult));
3669566063dSJacob Faibussowitsch       PetscCall(MatSetOption(aircraft.H,MAT_SYMMETRIC,PETSC_TRUE));
3679566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao,aircraft.H,aircraft.H,MatrixFreeObjHessian,(void*)&aircraft));
368c4762a1bSJed Brown     } else {
3699566063dSJacob Faibussowitsch       PetscCall(MatCreateDense(MPI_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,2*aircraft.nsteps,2*aircraft.nsteps,NULL,&(aircraft.H)));
3709566063dSJacob Faibussowitsch       PetscCall(MatSetOption(aircraft.H,MAT_SYMMETRIC,PETSC_TRUE));
3719566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao,aircraft.H,aircraft.H,FormObjHessian,(void *)&aircraft));
372c4762a1bSJed Brown     }
373c4762a1bSJed Brown   }
374c4762a1bSJed Brown 
375c4762a1bSJed Brown   /* Check for any TAO command line options */
3769566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao,&ksp));
377c4762a1bSJed Brown   if (ksp) {
3789566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp,&pc));
3799566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc,PCNONE));
380c4762a1bSJed Brown   }
3819566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
382c4762a1bSJed Brown 
3839566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
3849566063dSJacob Faibussowitsch   PetscCall(VecView(P,PETSC_VIEWER_STDOUT_WORLD));
385c4762a1bSJed Brown 
386c4762a1bSJed Brown   /* Free TAO data structures */
3879566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
388c4762a1bSJed Brown 
389c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
390c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
391c4762a1bSJed Brown      are no longer needed.
392c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3939566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&aircraft.ts));
3949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&aircraft.A));
3959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aircraft.U));
3969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aircraft.V));
3979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aircraft.W));
3989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&P));
3999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&PU));
4009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&PL));
4019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&aircraft.Jacp));
4029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&aircraft.DRDU));
4039566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&aircraft.DRDP));
4049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aircraft.Lambda[0]));
4059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aircraft.Mup[0]));
4069566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&P));
407c4762a1bSJed Brown   if (aircraft.eh) {
4089566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.Lambda2[0]));
4099566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.Mup2[0]));
4109566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.Dir));
4119566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.rhshp1[0]));
4129566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.rhshp2[0]));
4139566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.rhshp3[0]));
4149566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.rhshp4[0]));
4159566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.inthp1[0]));
4169566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.inthp2[0]));
4179566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.inthp3[0]));
4189566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.inthp4[0]));
4199566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&aircraft.H));
420c4762a1bSJed Brown   }
4219566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
422b122ec5aSJacob Faibussowitsch   return 0;
423c4762a1bSJed Brown }
424c4762a1bSJed Brown 
425c4762a1bSJed Brown /*
426c4762a1bSJed Brown    FormObjFunctionGradient - Evaluates the function and corresponding gradient.
427c4762a1bSJed Brown 
428c4762a1bSJed Brown    Input Parameters:
429c4762a1bSJed Brown    tao - the Tao context
430c4762a1bSJed Brown    P   - the input vector
431a82e8c82SStefano Zampini    ctx - optional aircraft-defined context, as set by TaoSetObjectiveAndGradient()
432c4762a1bSJed Brown 
433c4762a1bSJed Brown    Output Parameters:
434c4762a1bSJed Brown    f   - the newly evaluated function
435c4762a1bSJed Brown    G   - the newly evaluated gradient
436c4762a1bSJed Brown */
437c4762a1bSJed Brown PetscErrorCode FormObjFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
438c4762a1bSJed Brown {
439c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
440c4762a1bSJed Brown   TS                ts = actx->ts;
441c4762a1bSJed Brown   Vec               Q;
442c4762a1bSJed Brown   const PetscScalar *p,*q;
443c4762a1bSJed Brown   PetscScalar       *u,*v,*w;
444c4762a1bSJed Brown   PetscInt          i;
445c4762a1bSJed Brown 
446c4762a1bSJed Brown   PetscFunctionBeginUser;
4479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P,&p));
4489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->V,&v));
4499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->W,&w));
450c4762a1bSJed Brown   for (i=0; i<actx->nsteps; i++) {
451c4762a1bSJed Brown     v[i] = p[2*i];
452c4762a1bSJed Brown     w[i] = p[2*i+1];
453c4762a1bSJed Brown   }
4549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P,&p));
4559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->V,&v));
4569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->W,&w));
457c4762a1bSJed Brown 
4589566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts,0.0));
4599566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts,0));
4609566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
4619566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,actx->ftime/actx->nsteps));
462c4762a1bSJed Brown 
463c4762a1bSJed Brown   /* reinitialize system state */
4649566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->U,&u));
465c4762a1bSJed Brown   u[0] = 2.0;
466c4762a1bSJed Brown   u[1] = 0;
4679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->U,&u));
468c4762a1bSJed Brown 
469c4762a1bSJed Brown   /* reinitialize the integral value */
4709566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts,&Q));
4719566063dSJacob Faibussowitsch   PetscCall(VecSet(Q,0.0));
472c4762a1bSJed Brown 
4739566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,actx->U));
474c4762a1bSJed Brown 
475c4762a1bSJed Brown   /* Reset initial conditions for the adjoint integration */
4769566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Lambda[0],0.0));
4779566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Mup[0],0.0));
4789566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts,1,actx->Lambda,actx->Mup));
479c4762a1bSJed Brown 
4809566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
4819566063dSJacob Faibussowitsch   PetscCall(VecCopy(actx->Mup[0],G));
4829566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts,&Q));
4839566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Q,&q));
484c4762a1bSJed Brown   *f   = q[0];
4859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Q,&q));
486c4762a1bSJed Brown   PetscFunctionReturn(0);
487c4762a1bSJed Brown }
488c4762a1bSJed Brown 
489c4762a1bSJed Brown PetscErrorCode FormObjHessian(Tao tao,Vec P,Mat H,Mat Hpre,void *ctx)
490c4762a1bSJed Brown {
491c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
492c4762a1bSJed Brown   const PetscScalar *p;
493c4762a1bSJed Brown   PetscScalar       *harr,*v,*w,one = 1.0;
494c4762a1bSJed Brown   PetscInt          ind[1];
495c4762a1bSJed Brown   PetscInt          *cols,i;
496c4762a1bSJed Brown   Vec               Dir;
497c4762a1bSJed Brown 
498c4762a1bSJed Brown   PetscFunctionBeginUser;
499c4762a1bSJed Brown   /* set up control parameters */
5009566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P,&p));
5019566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->V,&v));
5029566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->W,&w));
503c4762a1bSJed Brown   for (i=0; i<actx->nsteps; i++) {
504c4762a1bSJed Brown     v[i] = p[2*i];
505c4762a1bSJed Brown     w[i] = p[2*i+1];
506c4762a1bSJed Brown   }
5079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P,&p));
5089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->V,&v));
5099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->W,&w));
510c4762a1bSJed Brown 
5119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*actx->nsteps,&harr));
5129566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2*actx->nsteps,&cols));
513c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) cols[i] = i;
5149566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(P,&Dir));
515c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) {
516c4762a1bSJed Brown     ind[0] = i;
5179566063dSJacob Faibussowitsch     PetscCall(VecSet(Dir,0.0));
5189566063dSJacob Faibussowitsch     PetscCall(VecSetValues(Dir,1,ind,&one,INSERT_VALUES));
5199566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(Dir));
5209566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(Dir));
5219566063dSJacob Faibussowitsch     PetscCall(ComputeObjHessianWithSOA(Dir,harr,actx));
5229566063dSJacob Faibussowitsch     PetscCall(MatSetValues(H,1,ind,2*actx->nsteps,cols,harr,INSERT_VALUES));
5239566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
5249566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
525c4762a1bSJed Brown     if (H != Hpre) {
5269566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(Hpre,MAT_FINAL_ASSEMBLY));
5279566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(Hpre,MAT_FINAL_ASSEMBLY));
528c4762a1bSJed Brown     }
529c4762a1bSJed Brown   }
5309566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
5319566063dSJacob Faibussowitsch   PetscCall(PetscFree(harr));
5329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Dir));
533c4762a1bSJed Brown   PetscFunctionReturn(0);
534c4762a1bSJed Brown }
535c4762a1bSJed Brown 
536c4762a1bSJed Brown PetscErrorCode MatrixFreeObjHessian(Tao tao, Vec P, Mat H, Mat Hpre, void *ctx)
537c4762a1bSJed Brown {
538c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
539c4762a1bSJed Brown   PetscScalar       *v,*w;
540c4762a1bSJed Brown   const PetscScalar *p;
541c4762a1bSJed Brown   PetscInt          i;
542c4762a1bSJed Brown 
543c4762a1bSJed Brown   PetscFunctionBegin;
5449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P,&p));
5459566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->V,&v));
5469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->W,&w));
547c4762a1bSJed Brown   for (i=0; i<actx->nsteps; i++) {
548c4762a1bSJed Brown     v[i] = p[2*i];
549c4762a1bSJed Brown     w[i] = p[2*i+1];
550c4762a1bSJed Brown   }
5519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P,&p));
5529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->V,&v));
5539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->W,&w));
554c4762a1bSJed Brown   PetscFunctionReturn(0);
555c4762a1bSJed Brown }
556c4762a1bSJed Brown 
557c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y)
558c4762a1bSJed Brown {
559c4762a1bSJed Brown   PetscScalar    *y;
560c4762a1bSJed Brown   void           *ptr;
561c4762a1bSJed Brown 
562c4762a1bSJed Brown   PetscFunctionBegin;
5639566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(H_shell,&ptr));
5649566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Y,&y));
5659566063dSJacob Faibussowitsch   PetscCall(ComputeObjHessianWithSOA(X,y,(Aircraft)ptr));
5669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Y,&y));
567c4762a1bSJed Brown   PetscFunctionReturn(0);
568c4762a1bSJed Brown }
569c4762a1bSJed Brown 
570c4762a1bSJed Brown PetscErrorCode ComputeObjHessianWithSOA(Vec Dir,PetscScalar arr[],Aircraft actx)
571c4762a1bSJed Brown {
572c4762a1bSJed Brown   TS                ts = actx->ts;
573c4762a1bSJed Brown   const PetscScalar *z_ptr;
574c4762a1bSJed Brown   PetscScalar       *u;
575c4762a1bSJed Brown   Vec               Q;
576c4762a1bSJed Brown   PetscInt          i;
577c4762a1bSJed Brown 
578c4762a1bSJed Brown   PetscFunctionBeginUser;
579c4762a1bSJed Brown   /* Reset TSAdjoint so that AdjointSetUp will be called again */
5809566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
581c4762a1bSJed Brown 
5829566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts,0.0));
5839566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts,0));
5849566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
5859566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,actx->ftime/actx->nsteps));
5869566063dSJacob Faibussowitsch   PetscCall(TSSetCostHessianProducts(actx->ts,1,actx->Lambda2,actx->Mup2,Dir));
587c4762a1bSJed Brown 
588c4762a1bSJed Brown   /* reinitialize system state */
5899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->U,&u));
590c4762a1bSJed Brown   u[0] = 2.0;
591c4762a1bSJed Brown   u[1] = 0;
5929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->U,&u));
593c4762a1bSJed Brown 
594c4762a1bSJed Brown   /* reinitialize the integral value */
5959566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts,&Q));
5969566063dSJacob Faibussowitsch   PetscCall(VecSet(Q,0.0));
597c4762a1bSJed Brown 
598c4762a1bSJed Brown   /* initialize tlm variable */
5999566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(actx->Jacp));
6009566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(actx->Jacp,MAT_FINAL_ASSEMBLY));
6019566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(actx->Jacp,MAT_FINAL_ASSEMBLY));
6029566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetForward(ts,actx->Jacp));
603c4762a1bSJed Brown 
6049566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,actx->U));
605c4762a1bSJed Brown 
606c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
6079566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Lambda[0],0.0));
6089566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Mup[0],0.0));
6099566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Lambda2[0],0.0));
6109566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Mup2[0],0.0));
6119566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts,1,actx->Lambda,actx->Mup));
612c4762a1bSJed Brown 
6139566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts,&Q));
614c4762a1bSJed Brown 
615c4762a1bSJed Brown   /* Reset initial conditions for the adjoint integration */
6169566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
617c4762a1bSJed Brown 
618c4762a1bSJed Brown   /* initial condition does not depend on p, so that lambda is not needed to assemble G */
6199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->Mup2[0],&z_ptr));
620c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) arr[i] = z_ptr[i];
6219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(actx->Mup2[0],&z_ptr));
622c4762a1bSJed Brown 
623c4762a1bSJed Brown   /* Disable second-order adjoint mode */
6249566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
6259566063dSJacob Faibussowitsch   PetscCall(TSAdjointResetForward(ts));
626c4762a1bSJed Brown   PetscFunctionReturn(0);
627c4762a1bSJed Brown }
628c4762a1bSJed Brown 
629c4762a1bSJed Brown /*TEST
630c4762a1bSJed Brown     build:
631c4762a1bSJed Brown       requires: !complex !single
632c4762a1bSJed Brown 
633c4762a1bSJed Brown     test:
634c4762a1bSJed Brown       args:  -ts_adapt_type none -ts_type rk -ts_rk_type 3 -viewer_binary_skip_info -tao_monitor -tao_gatol 1e-7
635c4762a1bSJed Brown 
636c4762a1bSJed Brown     test:
637c4762a1bSJed Brown       suffix: 2
638c4762a1bSJed 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
639c4762a1bSJed Brown 
640c4762a1bSJed Brown     test:
641c4762a1bSJed Brown       suffix: 3
642c4762a1bSJed 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
643c4762a1bSJed Brown TEST*/
644