xref: /petsc/src/ts/tutorials/optimal_control/ex1.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
325f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&step));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->V,&v));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->W,&w));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
37c4762a1bSJed Brown   f[0] = v[step]*PetscCosReal(w[step]);
38c4762a1bSJed Brown   f[1] = v[step]*PetscSinReal(w[step]);
395f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(actx->V,&v));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(actx->W,&w));
425f80ce2aSJacob 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;
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(A));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&step));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->V,&v));
585f80ce2aSJacob 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 
655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(actx->V,&v));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(actx->W,&w));
68c4762a1bSJed Brown 
697b0e2f17SHong Zhang   rowcol[0] = 2*step;
707b0e2f17SHong Zhang   rowcol[1] = 2*step+1;
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,2,rows,2,rowcol,&Jp[0][0],INSERT_VALUES));
72c4762a1bSJed Brown 
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
745f80ce2aSJacob 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;
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&step));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->V,&v));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->W,&w));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vl[0],&vl));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(VHV[0],0.0));
1125f80ce2aSJacob 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   }
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vl[0],&vl));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
1305f80ce2aSJacob 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;
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&step));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->V,&v));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->W,&w));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Vr,&vr));
1495f80ce2aSJacob 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   }
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Vr,&vr));
1615f80ce2aSJacob 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;
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
1925f80ce2aSJacob 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;
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(R,&r));
1975f80ce2aSJacob 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;
2105f80ce2aSJacob 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;
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(DRDU,2,row,1,col,&drdu[0][0],INSERT_VALUES));
2165f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY));
2185f80ce2aSJacob 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;
2255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(DRDP));
2265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY));
2275f80ce2aSJacob 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 
242c4762a1bSJed Brown   /* Initialize program */
243*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,NULL));
2445f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
2453c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /* Parameter settings */
248c4762a1bSJed Brown   aircraft.ftime = 1.;   /* time interval in hour */
249c4762a1bSJed Brown   aircraft.nsteps = 10; /* number of steps */
250c4762a1bSJed Brown   aircraft.lv = 2.0; /* leader speed in kmph */
251c4762a1bSJed Brown   aircraft.lw = PETSC_PI/4.; /* leader heading angle */
252c4762a1bSJed Brown 
2535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-ftime",&aircraft.ftime,NULL));
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nsteps",&aircraft.nsteps,NULL));
2555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-matrixfree",&aircraft.mf));
2565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-exacthessian",&aircraft.eh));
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
2595f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
2605f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOBQNLS));
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
2635f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&aircraft.A));
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(aircraft.A,PETSC_DECIDE,PETSC_DECIDE,2,2));
2655f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(aircraft.A));
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(aircraft.A));
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(aircraft.A,MAT_FINAL_ASSEMBLY));
2685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(aircraft.A,MAT_FINAL_ASSEMBLY));
2695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(aircraft.A,1));
2705f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(aircraft.A,-1));
271c4762a1bSJed Brown 
2725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&aircraft.Jacp));
2735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(aircraft.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,2*aircraft.nsteps));
2745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(aircraft.Jacp));
2755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(aircraft.Jacp));
2765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2*aircraft.nsteps,1,NULL,&aircraft.DRDP));
2775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(aircraft.DRDP));
2785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&aircraft.DRDU));
2795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(aircraft.DRDU));
280c4762a1bSJed Brown 
281c4762a1bSJed Brown   /* Create timestepping solver context */
2825f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&aircraft.ts));
2835f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(aircraft.ts,TSRK));
2845f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(aircraft.ts,NULL,RHSFunction,&aircraft));
2855f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(aircraft.ts,aircraft.A,aircraft.A,TSComputeRHSJacobianConstant,&aircraft));
2865f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobianP(aircraft.ts,aircraft.Jacp,RHSJacobianP,&aircraft));
2875f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(aircraft.ts,TS_EXACTFINALTIME_MATCHSTEP));
2885f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetEquationType(aircraft.ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   /* Set initial conditions */
2915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(aircraft.A,&aircraft.U,NULL));
2925f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(aircraft.ts,aircraft.U));
2935f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(aircraft.U,&u));
294c4762a1bSJed Brown   u[0] = 1.5;
295c4762a1bSJed Brown   u[1] = 0;
2965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(aircraft.U,&u));
2975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&aircraft.V));
2985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(aircraft.V,PETSC_DECIDE,aircraft.nsteps));
2995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetUp(aircraft.V));
3005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(aircraft.V,&aircraft.W));
3015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(aircraft.V,1.));
3025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(aircraft.W,PETSC_PI/4.));
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used */
3055f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSaveTrajectory(aircraft.ts));
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   /* Set sensitivity context */
3085f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreateQuadratureTS(aircraft.ts,PETSC_FALSE,&aircraft.quadts));
3095f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(aircraft.quadts,NULL,(TSRHSFunction)CostIntegrand,&aircraft));
3105f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(aircraft.quadts,aircraft.DRDU,aircraft.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&aircraft));
3115f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobianP(aircraft.quadts,aircraft.DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&aircraft));
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(aircraft.A,&aircraft.Lambda[0],NULL));
3135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(aircraft.Jacp,&aircraft.Mup[0],NULL));
314c4762a1bSJed Brown   if (aircraft.eh) {
3155f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.A,&aircraft.rhshp1[0],NULL));
3165f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.A,&aircraft.rhshp2[0],NULL));
3175f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.Jacp,&aircraft.rhshp3[0],NULL));
3185f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.Jacp,&aircraft.rhshp4[0],NULL));
3195f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.DRDU,&aircraft.inthp1[0],NULL));
3205f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.DRDU,&aircraft.inthp2[0],NULL));
3215f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.DRDP,&aircraft.inthp3[0],NULL));
3225f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.DRDP,&aircraft.inthp4[0],NULL));
3235f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.Jacp,&aircraft.Dir,NULL));
3245f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSHessianProduct(aircraft.ts,aircraft.rhshp1,RHSHessianProductUU,aircraft.rhshp2,RHSHessianProductUP,aircraft.rhshp3,RHSHessianProductPU,aircraft.rhshp4,RHSHessianProductPP,&aircraft));
3255f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSHessianProduct(aircraft.quadts,aircraft.inthp1,IntegrandHessianProductUU,aircraft.inthp2,IntegrandHessianProductUP,aircraft.inthp3,IntegrandHessianProductPU,aircraft.inthp4,IntegrandHessianProductPP,&aircraft));
3265f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.A,&aircraft.Lambda2[0],NULL));
3275f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(aircraft.Jacp,&aircraft.Mup2[0],NULL));
328c4762a1bSJed Brown   }
3295f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(aircraft.ts));
3305f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(aircraft.ts,aircraft.ftime));
3315f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(aircraft.ts,aircraft.ftime/aircraft.nsteps));
332c4762a1bSJed Brown 
333c4762a1bSJed Brown   /* Set initial solution guess */
3345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(aircraft.Jacp,&P,NULL));
3355f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(P,&p));
336c4762a1bSJed Brown   for (i=0; i<aircraft.nsteps; i++) {
337c4762a1bSJed Brown     p[2*i] = 2.0;
338c4762a1bSJed Brown     p[2*i+1] = PETSC_PI/2.0;
339c4762a1bSJed Brown   }
3405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(P,&p));
3415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(P,&PU));
3425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(P,&PL));
3435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(PU,&p));
344c4762a1bSJed Brown   for (i=0; i<aircraft.nsteps; i++) {
345c4762a1bSJed Brown     p[2*i] = 2.0;
346c4762a1bSJed Brown     p[2*i+1] = PETSC_PI;
347c4762a1bSJed Brown   }
3485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(PU,&p));
3495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(PL,&p));
350c4762a1bSJed Brown   for (i=0; i<aircraft.nsteps; i++) {
351c4762a1bSJed Brown     p[2*i] = 0.0;
352c4762a1bSJed Brown     p[2*i+1] = -PETSC_PI;
353c4762a1bSJed Brown   }
3545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(PL,&p));
355c4762a1bSJed Brown 
3565f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,P));
3575f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetVariableBounds(tao,PL,PU));
358c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
3595f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormObjFunctionGradient,(void *)&aircraft));
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   if (aircraft.eh) {
362c4762a1bSJed Brown     if (aircraft.mf) {
3635f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2*aircraft.nsteps,2*aircraft.nsteps,(void*)&aircraft,&aircraft.H));
3645f80ce2aSJacob Faibussowitsch       CHKERRQ(MatShellSetOperation(aircraft.H,MATOP_MULT,(void(*)(void))MyMatMult));
3655f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOption(aircraft.H,MAT_SYMMETRIC,PETSC_TRUE));
3665f80ce2aSJacob Faibussowitsch       CHKERRQ(TaoSetHessian(tao,aircraft.H,aircraft.H,MatrixFreeObjHessian,(void*)&aircraft));
367c4762a1bSJed Brown     } else {
3685f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateDense(MPI_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,2*aircraft.nsteps,2*aircraft.nsteps,NULL,&(aircraft.H)));
3695f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetOption(aircraft.H,MAT_SYMMETRIC,PETSC_TRUE));
3705f80ce2aSJacob Faibussowitsch       CHKERRQ(TaoSetHessian(tao,aircraft.H,aircraft.H,FormObjHessian,(void *)&aircraft));
371c4762a1bSJed Brown     }
372c4762a1bSJed Brown   }
373c4762a1bSJed Brown 
374c4762a1bSJed Brown   /* Check for any TAO command line options */
3755f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetKSP(tao,&ksp));
376c4762a1bSJed Brown   if (ksp) {
3775f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp,&pc));
3785f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetType(pc,PCNONE));
379c4762a1bSJed Brown   }
3805f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
381c4762a1bSJed Brown 
3825f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
3835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(P,PETSC_VIEWER_STDOUT_WORLD));
384c4762a1bSJed Brown 
385c4762a1bSJed Brown   /* Free TAO data structures */
3865f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
389c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
390c4762a1bSJed Brown      are no longer needed.
391c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3925f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&aircraft.ts));
3935f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&aircraft.A));
3945f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&aircraft.U));
3955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&aircraft.V));
3965f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&aircraft.W));
3975f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&P));
3985f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&PU));
3995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&PL));
4005f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&aircraft.Jacp));
4015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&aircraft.DRDU));
4025f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&aircraft.DRDP));
4035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&aircraft.Lambda[0]));
4045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&aircraft.Mup[0]));
4055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&P));
406c4762a1bSJed Brown   if (aircraft.eh) {
4075f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.Lambda2[0]));
4085f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.Mup2[0]));
4095f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.Dir));
4105f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.rhshp1[0]));
4115f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.rhshp2[0]));
4125f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.rhshp3[0]));
4135f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.rhshp4[0]));
4145f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.inthp1[0]));
4155f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.inthp2[0]));
4165f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.inthp3[0]));
4175f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&aircraft.inthp4[0]));
4185f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&aircraft.H));
419c4762a1bSJed Brown   }
420*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
421*b122ec5aSJacob Faibussowitsch   return 0;
422c4762a1bSJed Brown }
423c4762a1bSJed Brown 
424c4762a1bSJed Brown /*
425c4762a1bSJed Brown    FormObjFunctionGradient - Evaluates the function and corresponding gradient.
426c4762a1bSJed Brown 
427c4762a1bSJed Brown    Input Parameters:
428c4762a1bSJed Brown    tao - the Tao context
429c4762a1bSJed Brown    P   - the input vector
430a82e8c82SStefano Zampini    ctx - optional aircraft-defined context, as set by TaoSetObjectiveAndGradient()
431c4762a1bSJed Brown 
432c4762a1bSJed Brown    Output Parameters:
433c4762a1bSJed Brown    f   - the newly evaluated function
434c4762a1bSJed Brown    G   - the newly evaluated gradient
435c4762a1bSJed Brown */
436c4762a1bSJed Brown PetscErrorCode FormObjFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
437c4762a1bSJed Brown {
438c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
439c4762a1bSJed Brown   TS                ts = actx->ts;
440c4762a1bSJed Brown   Vec               Q;
441c4762a1bSJed Brown   const PetscScalar *p,*q;
442c4762a1bSJed Brown   PetscScalar       *u,*v,*w;
443c4762a1bSJed Brown   PetscInt          i;
444c4762a1bSJed Brown 
445c4762a1bSJed Brown   PetscFunctionBeginUser;
4465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,&p));
4475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->V,&v));
4485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->W,&w));
449c4762a1bSJed Brown   for (i=0; i<actx->nsteps; i++) {
450c4762a1bSJed Brown     v[i] = p[2*i];
451c4762a1bSJed Brown     w[i] = p[2*i+1];
452c4762a1bSJed Brown   }
4535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,&p));
4545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->V,&v));
4555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->W,&w));
456c4762a1bSJed Brown 
4575f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0.0));
4585f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ts,0));
4595f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
4605f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,actx->ftime/actx->nsteps));
461c4762a1bSJed Brown 
462c4762a1bSJed Brown   /* reinitialize system state */
4635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->U,&u));
464c4762a1bSJed Brown   u[0] = 2.0;
465c4762a1bSJed Brown   u[1] = 0;
4665f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->U,&u));
467c4762a1bSJed Brown 
468c4762a1bSJed Brown   /* reinitialize the integral value */
4695f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetCostIntegral(ts,&Q));
4705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(Q,0.0));
471c4762a1bSJed Brown 
4725f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,actx->U));
473c4762a1bSJed Brown 
474c4762a1bSJed Brown   /* Reset initial conditions for the adjoint integration */
4755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Lambda[0],0.0));
4765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Mup[0],0.0));
4775f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,1,actx->Lambda,actx->Mup));
478c4762a1bSJed Brown 
4795f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSolve(ts));
4805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(actx->Mup[0],G));
4815f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetCostIntegral(ts,&Q));
4825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Q,&q));
483c4762a1bSJed Brown   *f   = q[0];
4845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Q,&q));
485c4762a1bSJed Brown   PetscFunctionReturn(0);
486c4762a1bSJed Brown }
487c4762a1bSJed Brown 
488c4762a1bSJed Brown PetscErrorCode FormObjHessian(Tao tao,Vec P,Mat H,Mat Hpre,void *ctx)
489c4762a1bSJed Brown {
490c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
491c4762a1bSJed Brown   const PetscScalar *p;
492c4762a1bSJed Brown   PetscScalar       *harr,*v,*w,one = 1.0;
493c4762a1bSJed Brown   PetscInt          ind[1];
494c4762a1bSJed Brown   PetscInt          *cols,i;
495c4762a1bSJed Brown   Vec               Dir;
496c4762a1bSJed Brown 
497c4762a1bSJed Brown   PetscFunctionBeginUser;
498c4762a1bSJed Brown   /* set up control parameters */
4995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,&p));
5005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->V,&v));
5015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->W,&w));
502c4762a1bSJed Brown   for (i=0; i<actx->nsteps; i++) {
503c4762a1bSJed Brown     v[i] = p[2*i];
504c4762a1bSJed Brown     w[i] = p[2*i+1];
505c4762a1bSJed Brown   }
5065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,&p));
5075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->V,&v));
5085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->W,&w));
509c4762a1bSJed Brown 
5105f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(2*actx->nsteps,&harr));
5115f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(2*actx->nsteps,&cols));
512c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) cols[i] = i;
5135f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(P,&Dir));
514c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) {
515c4762a1bSJed Brown     ind[0] = i;
5165f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSet(Dir,0.0));
5175f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(Dir,1,ind,&one,INSERT_VALUES));
5185f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAssemblyBegin(Dir));
5195f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAssemblyEnd(Dir));
5205f80ce2aSJacob Faibussowitsch     CHKERRQ(ComputeObjHessianWithSOA(Dir,harr,actx));
5215f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(H,1,ind,2*actx->nsteps,cols,harr,INSERT_VALUES));
5225f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
5235f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
524c4762a1bSJed Brown     if (H != Hpre) {
5255f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyBegin(Hpre,MAT_FINAL_ASSEMBLY));
5265f80ce2aSJacob Faibussowitsch       CHKERRQ(MatAssemblyEnd(Hpre,MAT_FINAL_ASSEMBLY));
527c4762a1bSJed Brown     }
528c4762a1bSJed Brown   }
5295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(cols));
5305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(harr));
5315f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&Dir));
532c4762a1bSJed Brown   PetscFunctionReturn(0);
533c4762a1bSJed Brown }
534c4762a1bSJed Brown 
535c4762a1bSJed Brown PetscErrorCode MatrixFreeObjHessian(Tao tao, Vec P, Mat H, Mat Hpre, void *ctx)
536c4762a1bSJed Brown {
537c4762a1bSJed Brown   Aircraft          actx = (Aircraft)ctx;
538c4762a1bSJed Brown   PetscScalar       *v,*w;
539c4762a1bSJed Brown   const PetscScalar *p;
540c4762a1bSJed Brown   PetscInt          i;
541c4762a1bSJed Brown 
542c4762a1bSJed Brown   PetscFunctionBegin;
5435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,&p));
5445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->V,&v));
5455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->W,&w));
546c4762a1bSJed Brown   for (i=0; i<actx->nsteps; i++) {
547c4762a1bSJed Brown     v[i] = p[2*i];
548c4762a1bSJed Brown     w[i] = p[2*i+1];
549c4762a1bSJed Brown   }
5505f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,&p));
5515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->V,&v));
5525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->W,&w));
553c4762a1bSJed Brown   PetscFunctionReturn(0);
554c4762a1bSJed Brown }
555c4762a1bSJed Brown 
556c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y)
557c4762a1bSJed Brown {
558c4762a1bSJed Brown   PetscScalar    *y;
559c4762a1bSJed Brown   void           *ptr;
560c4762a1bSJed Brown 
561c4762a1bSJed Brown   PetscFunctionBegin;
5625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(H_shell,&ptr));
5635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(Y,&y));
5645f80ce2aSJacob Faibussowitsch   CHKERRQ(ComputeObjHessianWithSOA(X,y,(Aircraft)ptr));
5655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(Y,&y));
566c4762a1bSJed Brown   PetscFunctionReturn(0);
567c4762a1bSJed Brown }
568c4762a1bSJed Brown 
569c4762a1bSJed Brown PetscErrorCode ComputeObjHessianWithSOA(Vec Dir,PetscScalar arr[],Aircraft actx)
570c4762a1bSJed Brown {
571c4762a1bSJed Brown   TS                ts = actx->ts;
572c4762a1bSJed Brown   const PetscScalar *z_ptr;
573c4762a1bSJed Brown   PetscScalar       *u;
574c4762a1bSJed Brown   Vec               Q;
575c4762a1bSJed Brown   PetscInt          i;
576c4762a1bSJed Brown 
577c4762a1bSJed Brown   PetscFunctionBeginUser;
578c4762a1bSJed Brown   /* Reset TSAdjoint so that AdjointSetUp will be called again */
5795f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointReset(ts));
580c4762a1bSJed Brown 
5815f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0.0));
5825f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ts,0));
5835f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
5845f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,actx->ftime/actx->nsteps));
5855f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostHessianProducts(actx->ts,1,actx->Lambda2,actx->Mup2,Dir));
586c4762a1bSJed Brown 
587c4762a1bSJed Brown   /* reinitialize system state */
5885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(actx->U,&u));
589c4762a1bSJed Brown   u[0] = 2.0;
590c4762a1bSJed Brown   u[1] = 0;
5915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(actx->U,&u));
592c4762a1bSJed Brown 
593c4762a1bSJed Brown   /* reinitialize the integral value */
5945f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetCostIntegral(ts,&Q));
5955f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(Q,0.0));
596c4762a1bSJed Brown 
597c4762a1bSJed Brown   /* initialize tlm variable */
5985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(actx->Jacp));
5995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(actx->Jacp,MAT_FINAL_ASSEMBLY));
6005f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(actx->Jacp,MAT_FINAL_ASSEMBLY));
6015f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSetForward(ts,actx->Jacp));
602c4762a1bSJed Brown 
6035f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,actx->U));
604c4762a1bSJed Brown 
605c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
6065f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Lambda[0],0.0));
6075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Mup[0],0.0));
6085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Lambda2[0],0.0));
6095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(actx->Mup2[0],0.0));
6105f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,1,actx->Lambda,actx->Mup));
611c4762a1bSJed Brown 
6125f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetCostIntegral(ts,&Q));
613c4762a1bSJed Brown 
614c4762a1bSJed Brown   /* Reset initial conditions for the adjoint integration */
6155f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSolve(ts));
616c4762a1bSJed Brown 
617c4762a1bSJed Brown   /* initial condition does not depend on p, so that lambda is not needed to assemble G */
6185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(actx->Mup2[0],&z_ptr));
619c4762a1bSJed Brown   for (i=0; i<2*actx->nsteps; i++) arr[i] = z_ptr[i];
6205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(actx->Mup2[0],&z_ptr));
621c4762a1bSJed Brown 
622c4762a1bSJed Brown   /* Disable second-order adjoint mode */
6235f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointReset(ts));
6245f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointResetForward(ts));
625c4762a1bSJed Brown   PetscFunctionReturn(0);
626c4762a1bSJed Brown }
627c4762a1bSJed Brown 
628c4762a1bSJed Brown /*TEST
629c4762a1bSJed Brown     build:
630c4762a1bSJed Brown       requires: !complex !single
631c4762a1bSJed Brown 
632c4762a1bSJed Brown     test:
633c4762a1bSJed Brown       args:  -ts_adapt_type none -ts_type rk -ts_rk_type 3 -viewer_binary_skip_info -tao_monitor -tao_gatol 1e-7
634c4762a1bSJed Brown 
635c4762a1bSJed Brown     test:
636c4762a1bSJed Brown       suffix: 2
637c4762a1bSJed 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
638c4762a1bSJed Brown 
639c4762a1bSJed Brown     test:
640c4762a1bSJed Brown       suffix: 3
641c4762a1bSJed 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
642c4762a1bSJed Brown TEST*/
643