xref: /petsc/src/ts/tutorials/optimal_control/ex1.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&step));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->V,&v));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->W,&w));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
37c4762a1bSJed Brown   f[0] = v[step]*PetscCosReal(w[step]);
38c4762a1bSJed Brown   f[1] = v[step]*PetscSinReal(w[step]);
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(actx->V,&v));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(actx->W,&w));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(A));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&step));
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->V,&v));
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(actx->V,&v));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(actx->W,&w));
68c4762a1bSJed Brown 
697b0e2f17SHong Zhang   rowcol[0] = 2*step;
707b0e2f17SHong Zhang   rowcol[1] = 2*step+1;
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,2,rows,2,rowcol,&Jp[0][0],INSERT_VALUES));
72c4762a1bSJed Brown 
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&step));
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->V,&v));
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->W,&w));
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
110*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
111*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(VHV[0],0.0));
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
130*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&step));
145*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->V,&v));
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->W,&w));
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
191*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
196*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(R,&r));
197*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
210*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
215*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(DRDU,2,row,1,col,&drdu[0][0],INSERT_VALUES));
216*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
217*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY));
218*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
225*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(DRDP));
226*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY));
227*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   PetscErrorCode     ierr;
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /* Initialize program */
244c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,NULL);if (ierr) return ierr;
245*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(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 
254*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-ftime",&aircraft.ftime,NULL));
255*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nsteps",&aircraft.nsteps,NULL));
256*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-matrixfree",&aircraft.mf));
257*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-exacthessian",&aircraft.eh));
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
260*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
261*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOBQNLS));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
264*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&aircraft.A));
265*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(aircraft.A,PETSC_DECIDE,PETSC_DECIDE,2,2));
266*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(aircraft.A));
267*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(aircraft.A));
268*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(aircraft.A,MAT_FINAL_ASSEMBLY));
269*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(aircraft.A,MAT_FINAL_ASSEMBLY));
270*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(aircraft.A,1));
271*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(aircraft.A,-1));
272c4762a1bSJed Brown 
273*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&aircraft.Jacp));
274*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(aircraft.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,2*aircraft.nsteps));
275*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(aircraft.Jacp));
276*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(aircraft.Jacp));
277*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2*aircraft.nsteps,1,NULL,&aircraft.DRDP));
278*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(aircraft.DRDP));
279*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&aircraft.DRDU));
280*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(aircraft.DRDU));
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /* Create timestepping solver context */
283*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&aircraft.ts));
284*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(aircraft.ts,TSRK));
285*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(aircraft.ts,NULL,RHSFunction,&aircraft));
286*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(aircraft.ts,aircraft.A,aircraft.A,TSComputeRHSJacobianConstant,&aircraft));
287*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobianP(aircraft.ts,aircraft.Jacp,RHSJacobianP,&aircraft));
288*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(aircraft.ts,TS_EXACTFINALTIME_MATCHSTEP));
289*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 */
292*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(aircraft.A,&aircraft.U,NULL));
293*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(aircraft.ts,aircraft.U));
294*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(aircraft.U,&u));
295c4762a1bSJed Brown   u[0] = 1.5;
296c4762a1bSJed Brown   u[1] = 0;
297*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(aircraft.U,&u));
298*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&aircraft.V));
299*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(aircraft.V,PETSC_DECIDE,aircraft.nsteps));
300*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetUp(aircraft.V));
301*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(aircraft.V,&aircraft.W));
302*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(aircraft.V,1.));
303*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(aircraft.W,PETSC_PI/4.));
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used */
306*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSaveTrajectory(aircraft.ts));
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   /* Set sensitivity context */
309*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreateQuadratureTS(aircraft.ts,PETSC_FALSE,&aircraft.quadts));
310*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(aircraft.quadts,NULL,(TSRHSFunction)CostIntegrand,&aircraft));
311*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(aircraft.quadts,aircraft.DRDU,aircraft.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&aircraft));
312*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobianP(aircraft.quadts,aircraft.DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&aircraft));
313*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(aircraft.A,&aircraft.Lambda[0],NULL));
314*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(aircraft.Jacp,&aircraft.Mup[0],NULL));
315c4762a1bSJed Brown   if (aircraft.eh) {
316*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.A,&aircraft.rhshp1[0],NULL));
317*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.A,&aircraft.rhshp2[0],NULL));
318*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.Jacp,&aircraft.rhshp3[0],NULL));
319*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.Jacp,&aircraft.rhshp4[0],NULL));
320*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.DRDU,&aircraft.inthp1[0],NULL));
321*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.DRDU,&aircraft.inthp2[0],NULL));
322*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.DRDP,&aircraft.inthp3[0],NULL));
323*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.DRDP,&aircraft.inthp4[0],NULL));
324*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.Jacp,&aircraft.Dir,NULL));
325*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSHessianProduct(aircraft.ts,aircraft.rhshp1,RHSHessianProductUU,aircraft.rhshp2,RHSHessianProductUP,aircraft.rhshp3,RHSHessianProductPU,aircraft.rhshp4,RHSHessianProductPP,&aircraft));
326*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSHessianProduct(aircraft.quadts,aircraft.inthp1,IntegrandHessianProductUU,aircraft.inthp2,IntegrandHessianProductUP,aircraft.inthp3,IntegrandHessianProductPU,aircraft.inthp4,IntegrandHessianProductPP,&aircraft));
327*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.A,&aircraft.Lambda2[0],NULL));
328*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.Jacp,&aircraft.Mup2[0],NULL));
329c4762a1bSJed Brown   }
330*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(aircraft.ts));
331*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(aircraft.ts,aircraft.ftime));
332*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(aircraft.ts,aircraft.ftime/aircraft.nsteps));
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   /* Set initial solution guess */
335*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(aircraft.Jacp,&P,NULL));
336*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
341*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(P,&p));
342*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(P,&PU));
343*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(P,&PL));
344*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
349*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(PU,&p));
350*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
355*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(PL,&p));
356c4762a1bSJed Brown 
357*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,P));
358*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetVariableBounds(tao,PL,PU));
359c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
360*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormObjFunctionGradient,(void *)&aircraft));
361c4762a1bSJed Brown 
362c4762a1bSJed Brown   if (aircraft.eh) {
363c4762a1bSJed Brown     if (aircraft.mf) {
364*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2*aircraft.nsteps,2*aircraft.nsteps,(void*)&aircraft,&aircraft.H));
365*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatShellSetOperation(aircraft.H,MATOP_MULT,(void(*)(void))MyMatMult));
366*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOption(aircraft.H,MAT_SYMMETRIC,PETSC_TRUE));
367*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TaoSetHessian(tao,aircraft.H,aircraft.H,MatrixFreeObjHessian,(void*)&aircraft));
368c4762a1bSJed Brown     } else {
369*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateDense(MPI_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,2*aircraft.nsteps,2*aircraft.nsteps,NULL,&(aircraft.H)));
370*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOption(aircraft.H,MAT_SYMMETRIC,PETSC_TRUE));
371*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TaoSetHessian(tao,aircraft.H,aircraft.H,FormObjHessian,(void *)&aircraft));
372c4762a1bSJed Brown     }
373c4762a1bSJed Brown   }
374c4762a1bSJed Brown 
375c4762a1bSJed Brown   /* Check for any TAO command line options */
376*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetKSP(tao,&ksp));
377c4762a1bSJed Brown   if (ksp) {
378*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp,&pc));
379*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetType(pc,PCNONE));
380c4762a1bSJed Brown   }
381*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
382c4762a1bSJed Brown 
383*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
384*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(P,PETSC_VIEWER_STDOUT_WORLD));
385c4762a1bSJed Brown 
386c4762a1bSJed Brown   /* Free TAO data structures */
387*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
393*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&aircraft.ts));
394*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&aircraft.A));
395*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&aircraft.U));
396*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&aircraft.V));
397*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&aircraft.W));
398*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&P));
399*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&PU));
400*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&PL));
401*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&aircraft.Jacp));
402*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&aircraft.DRDU));
403*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&aircraft.DRDP));
404*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&aircraft.Lambda[0]));
405*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&aircraft.Mup[0]));
406*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&P));
407c4762a1bSJed Brown   if (aircraft.eh) {
408*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.Lambda2[0]));
409*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.Mup2[0]));
410*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.Dir));
411*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.rhshp1[0]));
412*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.rhshp2[0]));
413*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.rhshp3[0]));
414*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.rhshp4[0]));
415*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.inthp1[0]));
416*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.inthp2[0]));
417*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.inthp3[0]));
418*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.inthp4[0]));
419*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&aircraft.H));
420c4762a1bSJed Brown   }
421c4762a1bSJed Brown   ierr = PetscFinalize();
422c4762a1bSJed Brown   return ierr;
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;
447*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,&p));
448*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->V,&v));
449*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
454*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,&p));
455*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->V,&v));
456*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->W,&w));
457c4762a1bSJed Brown 
458*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0.0));
459*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ts,0));
460*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
461*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,actx->ftime/actx->nsteps));
462c4762a1bSJed Brown 
463c4762a1bSJed Brown   /* reinitialize system state */
464*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->U,&u));
465c4762a1bSJed Brown   u[0] = 2.0;
466c4762a1bSJed Brown   u[1] = 0;
467*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->U,&u));
468c4762a1bSJed Brown 
469c4762a1bSJed Brown   /* reinitialize the integral value */
470*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetCostIntegral(ts,&Q));
471*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(Q,0.0));
472c4762a1bSJed Brown 
473*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,actx->U));
474c4762a1bSJed Brown 
475c4762a1bSJed Brown   /* Reset initial conditions for the adjoint integration */
476*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Lambda[0],0.0));
477*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Mup[0],0.0));
478*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,1,actx->Lambda,actx->Mup));
479c4762a1bSJed Brown 
480*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSolve(ts));
481*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(actx->Mup[0],G));
482*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetCostIntegral(ts,&Q));
483*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Q,&q));
484c4762a1bSJed Brown   *f   = q[0];
485*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 */
500*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,&p));
501*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->V,&v));
502*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
507*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,&p));
508*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->V,&v));
509*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->W,&w));
510c4762a1bSJed Brown 
511*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(2*actx->nsteps,&harr));
512*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(2*actx->nsteps,&cols));
513c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) cols[i] = i;
514*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(P,&Dir));
515c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) {
516c4762a1bSJed Brown     ind[0] = i;
517*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(Dir,0.0));
518*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(Dir,1,ind,&one,INSERT_VALUES));
519*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAssemblyBegin(Dir));
520*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAssemblyEnd(Dir));
521*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ComputeObjHessianWithSOA(Dir,harr,actx));
522*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(H,1,ind,2*actx->nsteps,cols,harr,INSERT_VALUES));
523*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
524*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
525c4762a1bSJed Brown     if (H != Hpre) {
526*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyBegin(Hpre,MAT_FINAL_ASSEMBLY));
527*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyEnd(Hpre,MAT_FINAL_ASSEMBLY));
528c4762a1bSJed Brown     }
529c4762a1bSJed Brown   }
530*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(cols));
531*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(harr));
532*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
544*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,&p));
545*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->V,&v));
546*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   }
551*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,&p));
552*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->V,&v));
553*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
563*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(H_shell,&ptr));
564*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(Y,&y));
565*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ComputeObjHessianWithSOA(X,y,(Aircraft)ptr));
566*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 */
580*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointReset(ts));
581c4762a1bSJed Brown 
582*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0.0));
583*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ts,0));
584*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
585*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,actx->ftime/actx->nsteps));
586*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostHessianProducts(actx->ts,1,actx->Lambda2,actx->Mup2,Dir));
587c4762a1bSJed Brown 
588c4762a1bSJed Brown   /* reinitialize system state */
589*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->U,&u));
590c4762a1bSJed Brown   u[0] = 2.0;
591c4762a1bSJed Brown   u[1] = 0;
592*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->U,&u));
593c4762a1bSJed Brown 
594c4762a1bSJed Brown   /* reinitialize the integral value */
595*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetCostIntegral(ts,&Q));
596*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(Q,0.0));
597c4762a1bSJed Brown 
598c4762a1bSJed Brown   /* initialize tlm variable */
599*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(actx->Jacp));
600*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(actx->Jacp,MAT_FINAL_ASSEMBLY));
601*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(actx->Jacp,MAT_FINAL_ASSEMBLY));
602*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSetForward(ts,actx->Jacp));
603c4762a1bSJed Brown 
604*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,actx->U));
605c4762a1bSJed Brown 
606c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
607*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Lambda[0],0.0));
608*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Mup[0],0.0));
609*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Lambda2[0],0.0));
610*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Mup2[0],0.0));
611*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,1,actx->Lambda,actx->Mup));
612c4762a1bSJed Brown 
613*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetCostIntegral(ts,&Q));
614c4762a1bSJed Brown 
615c4762a1bSJed Brown   /* Reset initial conditions for the adjoint integration */
616*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSolve(ts));
617c4762a1bSJed Brown 
618c4762a1bSJed Brown   /* initial condition does not depend on p, so that lambda is not needed to assemble G */
619*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->Mup2[0],&z_ptr));
620c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) arr[i] = z_ptr[i];
621*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(actx->Mup2[0],&z_ptr));
622c4762a1bSJed Brown 
623c4762a1bSJed Brown   /* Disable second-order adjoint mode */
624*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointReset(ts));
625*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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