xref: /petsc/src/ts/tutorials/optimal_control/ex1.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
24*9371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx) {
25c4762a1bSJed Brown   Aircraft           actx = (Aircraft)ctx;
26c4762a1bSJed Brown   const PetscScalar *u, *v, *w;
27c4762a1bSJed Brown   PetscScalar       *f;
28c4762a1bSJed Brown   PetscInt           step;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   PetscFunctionBeginUser;
319566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &step));
329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->V, &v));
349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->W, &w));
359566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
36c4762a1bSJed Brown   f[0] = v[step] * PetscCosReal(w[step]);
37c4762a1bSJed Brown   f[1] = v[step] * PetscSinReal(w[step]);
389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(actx->V, &v));
409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(actx->W, &w));
419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
42c4762a1bSJed Brown   PetscFunctionReturn(0);
43c4762a1bSJed Brown }
44c4762a1bSJed Brown 
45*9371c9d4SSatish Balay static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec U, Mat A, void *ctx) {
46c4762a1bSJed Brown   Aircraft           actx = (Aircraft)ctx;
47c4762a1bSJed Brown   const PetscScalar *u, *v, *w;
48c4762a1bSJed Brown   PetscInt           step, rows[2] = {0, 1}, rowcol[2];
49c4762a1bSJed Brown   PetscScalar        Jp[2][2];
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   PetscFunctionBeginUser;
529566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(A));
539566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &step));
549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->V, &v));
569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->W, &w));
57c4762a1bSJed Brown 
587b0e2f17SHong Zhang   Jp[0][0] = PetscCosReal(w[step]);
597b0e2f17SHong Zhang   Jp[0][1] = -v[step] * PetscSinReal(w[step]);
607b0e2f17SHong Zhang   Jp[1][0] = PetscSinReal(w[step]);
617b0e2f17SHong Zhang   Jp[1][1] = v[step] * PetscCosReal(w[step]);
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(actx->V, &v));
659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(actx->W, &w));
66c4762a1bSJed Brown 
677b0e2f17SHong Zhang   rowcol[0] = 2 * step;
687b0e2f17SHong Zhang   rowcol[1] = 2 * step + 1;
699566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, rows, 2, rowcol, &Jp[0][0], INSERT_VALUES));
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
729566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
73c4762a1bSJed Brown   PetscFunctionReturn(0);
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76*9371c9d4SSatish Balay static PetscErrorCode RHSHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
77c4762a1bSJed Brown   PetscFunctionBeginUser;
78c4762a1bSJed Brown   PetscFunctionReturn(0);
79c4762a1bSJed Brown }
80c4762a1bSJed Brown 
81*9371c9d4SSatish Balay static PetscErrorCode RHSHessianProductUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
82c4762a1bSJed Brown   PetscFunctionBeginUser;
83c4762a1bSJed Brown   PetscFunctionReturn(0);
84c4762a1bSJed Brown }
85c4762a1bSJed Brown 
86*9371c9d4SSatish Balay static PetscErrorCode RHSHessianProductPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
87c4762a1bSJed Brown   PetscFunctionBeginUser;
88c4762a1bSJed Brown   PetscFunctionReturn(0);
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91*9371c9d4SSatish Balay static PetscErrorCode RHSHessianProductPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
92c4762a1bSJed Brown   Aircraft           actx = (Aircraft)ctx;
93c4762a1bSJed Brown   const PetscScalar *v, *w, *vl, *vr, *u;
94c4762a1bSJed Brown   PetscScalar       *vhv;
95c4762a1bSJed Brown   PetscScalar        dJpdP[2][2][2] = {{{0}}};
96c4762a1bSJed Brown   PetscInt           step, i, j, k;
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   PetscFunctionBeginUser;
999566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &step));
1009566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1019566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->V, &v));
1029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->W, &w));
1039566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vl[0], &vl));
1049566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
1059566063dSJacob Faibussowitsch   PetscCall(VecSet(VHV[0], 0.0));
1069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
107c4762a1bSJed Brown 
1087b0e2f17SHong Zhang   dJpdP[0][0][1] = -PetscSinReal(w[step]);
1097b0e2f17SHong Zhang   dJpdP[0][1][0] = -PetscSinReal(w[step]);
1107b0e2f17SHong Zhang   dJpdP[0][1][1] = -v[step] * PetscCosReal(w[step]);
1117b0e2f17SHong Zhang   dJpdP[1][0][1] = PetscCosReal(w[step]);
1127b0e2f17SHong Zhang   dJpdP[1][1][0] = PetscCosReal(w[step]);
1137b0e2f17SHong Zhang   dJpdP[1][1][1] = -v[step] * PetscSinReal(w[step]);
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
1167b0e2f17SHong Zhang     vhv[2 * step + j] = 0;
117c4762a1bSJed Brown     for (k = 0; k < 2; k++)
118*9371c9d4SSatish Balay       for (i = 0; i < 2; i++) vhv[2 * step + j] += vl[i] * dJpdP[i][j][k] * vr[2 * step + k];
119c4762a1bSJed Brown   }
1209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vl[0], &vl));
1229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
124c4762a1bSJed Brown   PetscFunctionReturn(0);
125c4762a1bSJed Brown }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown /* Vl in NULL,updates to VHV must be added */
128*9371c9d4SSatish Balay static PetscErrorCode IntegrandHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
129c4762a1bSJed Brown   Aircraft           actx = (Aircraft)ctx;
130c4762a1bSJed Brown   const PetscScalar *v, *w, *vr, *u;
131c4762a1bSJed Brown   PetscScalar       *vhv;
132c4762a1bSJed Brown   PetscScalar        dRudU[2][2] = {{0}};
133c4762a1bSJed Brown   PetscInt           step, j, k;
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   PetscFunctionBeginUser;
1369566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &step));
1379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->V, &v));
1399566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->W, &w));
1409566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Vr, &vr));
1419566063dSJacob Faibussowitsch   PetscCall(VecGetArray(VHV[0], &vhv));
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   dRudU[0][0] = 2.0;
144c4762a1bSJed Brown   dRudU[1][1] = 2.0;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   for (j = 0; j < 2; j++) {
147c4762a1bSJed Brown     vhv[j] = 0;
148*9371c9d4SSatish Balay     for (k = 0; k < 2; k++) vhv[j] += dRudU[j][k] * vr[k];
149c4762a1bSJed Brown   }
1509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Vr, &vr));
1529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(VHV[0], &vhv));
153c4762a1bSJed Brown   PetscFunctionReturn(0);
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
156*9371c9d4SSatish Balay static PetscErrorCode IntegrandHessianProductUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
157c4762a1bSJed Brown   PetscFunctionBeginUser;
158c4762a1bSJed Brown   PetscFunctionReturn(0);
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161*9371c9d4SSatish Balay static PetscErrorCode IntegrandHessianProductPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
162c4762a1bSJed Brown   PetscFunctionBeginUser;
163c4762a1bSJed Brown   PetscFunctionReturn(0);
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
166*9371c9d4SSatish Balay static PetscErrorCode IntegrandHessianProductPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) {
167c4762a1bSJed Brown   PetscFunctionBeginUser;
168c4762a1bSJed Brown   PetscFunctionReturn(0);
169c4762a1bSJed Brown }
170c4762a1bSJed Brown 
171*9371c9d4SSatish Balay static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, void *ctx) {
172c4762a1bSJed Brown   Aircraft           actx = (Aircraft)ctx;
173c4762a1bSJed Brown   PetscScalar       *r;
174c4762a1bSJed Brown   PetscReal          dx, dy;
175c4762a1bSJed Brown   const PetscScalar *u;
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   PetscFunctionBegin;
1789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
180c4762a1bSJed Brown   dx   = u[0] - actx->lv * t * PetscCosReal(actx->lw);
181c4762a1bSJed Brown   dy   = u[1] - actx->lv * t * PetscSinReal(actx->lw);
182c4762a1bSJed Brown   r[0] = dx * dx + dy * dy;
1839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
1849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
185c4762a1bSJed Brown   PetscFunctionReturn(0);
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188*9371c9d4SSatish Balay static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, void *ctx) {
189c4762a1bSJed Brown   Aircraft           actx = (Aircraft)ctx;
190c4762a1bSJed Brown   PetscScalar        drdu[2][1];
191c4762a1bSJed Brown   const PetscScalar *u;
192c4762a1bSJed Brown   PetscReal          dx, dy;
193c4762a1bSJed Brown   PetscInt           row[] = {0, 1}, col[] = {0};
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   PetscFunctionBegin;
1969566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
197c4762a1bSJed Brown   dx         = u[0] - actx->lv * t * PetscCosReal(actx->lw);
198c4762a1bSJed Brown   dy         = u[1] - actx->lv * t * PetscSinReal(actx->lw);
199c4762a1bSJed Brown   drdu[0][0] = 2. * dx;
200c4762a1bSJed Brown   drdu[1][0] = 2. * dy;
2019566063dSJacob Faibussowitsch   PetscCall(MatSetValues(DRDU, 2, row, 1, col, &drdu[0][0], INSERT_VALUES));
2029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
2039566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY));
2049566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY));
205c4762a1bSJed Brown   PetscFunctionReturn(0);
206c4762a1bSJed Brown }
207c4762a1bSJed Brown 
208*9371c9d4SSatish Balay static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, void *ctx) {
209c4762a1bSJed Brown   PetscFunctionBegin;
2109566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(DRDP));
2119566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY));
2129566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY));
213c4762a1bSJed Brown   PetscFunctionReturn(0);
214c4762a1bSJed Brown }
215c4762a1bSJed Brown 
216*9371c9d4SSatish Balay int main(int argc, char **argv) {
217c4762a1bSJed Brown   Vec                P, PL, PU;
218c4762a1bSJed Brown   struct _n_aircraft aircraft;
219c4762a1bSJed Brown   PetscMPIInt        size;
220c4762a1bSJed Brown   Tao                tao;
221c4762a1bSJed Brown   KSP                ksp;
222c4762a1bSJed Brown   PC                 pc;
223c4762a1bSJed Brown   PetscScalar       *u, *p;
224c4762a1bSJed Brown   PetscInt           i;
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /* Initialize program */
227327415f7SBarry Smith   PetscFunctionBeginUser;
2289566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, NULL));
2299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2303c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* Parameter settings */
233c4762a1bSJed Brown   aircraft.ftime  = 1.;            /* time interval in hour */
234c4762a1bSJed Brown   aircraft.nsteps = 10;            /* number of steps */
235c4762a1bSJed Brown   aircraft.lv     = 2.0;           /* leader speed in kmph */
236c4762a1bSJed Brown   aircraft.lw     = PETSC_PI / 4.; /* leader heading angle */
237c4762a1bSJed Brown 
2389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-ftime", &aircraft.ftime, NULL));
2399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nsteps", &aircraft.nsteps, NULL));
2409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-matrixfree", &aircraft.mf));
2419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-exacthessian", &aircraft.eh));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
2449566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
2459566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOBQNLS));
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /* Create necessary matrix and vectors, solve same ODE on every process */
2489566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &aircraft.A));
2499566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(aircraft.A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
2509566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(aircraft.A));
2519566063dSJacob Faibussowitsch   PetscCall(MatSetUp(aircraft.A));
2529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(aircraft.A, MAT_FINAL_ASSEMBLY));
2539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(aircraft.A, MAT_FINAL_ASSEMBLY));
2549566063dSJacob Faibussowitsch   PetscCall(MatShift(aircraft.A, 1));
2559566063dSJacob Faibussowitsch   PetscCall(MatShift(aircraft.A, -1));
256c4762a1bSJed Brown 
2579566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &aircraft.Jacp));
2589566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(aircraft.Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 2 * aircraft.nsteps));
2599566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(aircraft.Jacp));
2609566063dSJacob Faibussowitsch   PetscCall(MatSetUp(aircraft.Jacp));
2619566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2 * aircraft.nsteps, 1, NULL, &aircraft.DRDP));
2629566063dSJacob Faibussowitsch   PetscCall(MatSetUp(aircraft.DRDP));
2639566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, &aircraft.DRDU));
2649566063dSJacob Faibussowitsch   PetscCall(MatSetUp(aircraft.DRDU));
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   /* Create timestepping solver context */
2679566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &aircraft.ts));
2689566063dSJacob Faibussowitsch   PetscCall(TSSetType(aircraft.ts, TSRK));
2699566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(aircraft.ts, NULL, RHSFunction, &aircraft));
2709566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(aircraft.ts, aircraft.A, aircraft.A, TSComputeRHSJacobianConstant, &aircraft));
2719566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(aircraft.ts, aircraft.Jacp, RHSJacobianP, &aircraft));
2729566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(aircraft.ts, TS_EXACTFINALTIME_MATCHSTEP));
2739566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(aircraft.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   /* Set initial conditions */
2769566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(aircraft.A, &aircraft.U, NULL));
2779566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(aircraft.ts, aircraft.U));
2789566063dSJacob Faibussowitsch   PetscCall(VecGetArray(aircraft.U, &u));
279c4762a1bSJed Brown   u[0] = 1.5;
280c4762a1bSJed Brown   u[1] = 0;
2819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(aircraft.U, &u));
2829566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &aircraft.V));
2839566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(aircraft.V, PETSC_DECIDE, aircraft.nsteps));
2849566063dSJacob Faibussowitsch   PetscCall(VecSetUp(aircraft.V));
2859566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(aircraft.V, &aircraft.W));
2869566063dSJacob Faibussowitsch   PetscCall(VecSet(aircraft.V, 1.));
2879566063dSJacob Faibussowitsch   PetscCall(VecSet(aircraft.W, PETSC_PI / 4.));
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   /* Save trajectory of solution so that TSAdjointSolve() may be used */
2909566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(aircraft.ts));
291c4762a1bSJed Brown 
292c4762a1bSJed Brown   /* Set sensitivity context */
2939566063dSJacob Faibussowitsch   PetscCall(TSCreateQuadratureTS(aircraft.ts, PETSC_FALSE, &aircraft.quadts));
2949566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(aircraft.quadts, NULL, (TSRHSFunction)CostIntegrand, &aircraft));
2959566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(aircraft.quadts, aircraft.DRDU, aircraft.DRDU, (TSRHSJacobian)DRDUJacobianTranspose, &aircraft));
2969566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(aircraft.quadts, aircraft.DRDP, (TSRHSJacobianP)DRDPJacobianTranspose, &aircraft));
2979566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(aircraft.A, &aircraft.Lambda[0], NULL));
2989566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(aircraft.Jacp, &aircraft.Mup[0], NULL));
299c4762a1bSJed Brown   if (aircraft.eh) {
3009566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.A, &aircraft.rhshp1[0], NULL));
3019566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.A, &aircraft.rhshp2[0], NULL));
3029566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.Jacp, &aircraft.rhshp3[0], NULL));
3039566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.Jacp, &aircraft.rhshp4[0], NULL));
3049566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.DRDU, &aircraft.inthp1[0], NULL));
3059566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.DRDU, &aircraft.inthp2[0], NULL));
3069566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.DRDP, &aircraft.inthp3[0], NULL));
3079566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.DRDP, &aircraft.inthp4[0], NULL));
3089566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.Jacp, &aircraft.Dir, NULL));
3099566063dSJacob Faibussowitsch     PetscCall(TSSetRHSHessianProduct(aircraft.ts, aircraft.rhshp1, RHSHessianProductUU, aircraft.rhshp2, RHSHessianProductUP, aircraft.rhshp3, RHSHessianProductPU, aircraft.rhshp4, RHSHessianProductPP, &aircraft));
3109566063dSJacob Faibussowitsch     PetscCall(TSSetRHSHessianProduct(aircraft.quadts, aircraft.inthp1, IntegrandHessianProductUU, aircraft.inthp2, IntegrandHessianProductUP, aircraft.inthp3, IntegrandHessianProductPU, aircraft.inthp4, IntegrandHessianProductPP, &aircraft));
3119566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.A, &aircraft.Lambda2[0], NULL));
3129566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(aircraft.Jacp, &aircraft.Mup2[0], NULL));
313c4762a1bSJed Brown   }
3149566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(aircraft.ts));
3159566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(aircraft.ts, aircraft.ftime));
3169566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(aircraft.ts, aircraft.ftime / aircraft.nsteps));
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   /* Set initial solution guess */
3199566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(aircraft.Jacp, &P, NULL));
3209566063dSJacob Faibussowitsch   PetscCall(VecGetArray(P, &p));
321c4762a1bSJed Brown   for (i = 0; i < aircraft.nsteps; i++) {
322c4762a1bSJed Brown     p[2 * i]     = 2.0;
323c4762a1bSJed Brown     p[2 * i + 1] = PETSC_PI / 2.0;
324c4762a1bSJed Brown   }
3259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(P, &p));
3269566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(P, &PU));
3279566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(P, &PL));
3289566063dSJacob Faibussowitsch   PetscCall(VecGetArray(PU, &p));
329c4762a1bSJed Brown   for (i = 0; i < aircraft.nsteps; i++) {
330c4762a1bSJed Brown     p[2 * i]     = 2.0;
331c4762a1bSJed Brown     p[2 * i + 1] = PETSC_PI;
332c4762a1bSJed Brown   }
3339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(PU, &p));
3349566063dSJacob Faibussowitsch   PetscCall(VecGetArray(PL, &p));
335c4762a1bSJed Brown   for (i = 0; i < aircraft.nsteps; i++) {
336c4762a1bSJed Brown     p[2 * i]     = 0.0;
337c4762a1bSJed Brown     p[2 * i + 1] = -PETSC_PI;
338c4762a1bSJed Brown   }
3399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(PL, &p));
340c4762a1bSJed Brown 
3419566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, P));
3429566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao, PL, PU));
343c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
3449566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormObjFunctionGradient, (void *)&aircraft));
345c4762a1bSJed Brown 
346c4762a1bSJed Brown   if (aircraft.eh) {
347c4762a1bSJed Brown     if (aircraft.mf) {
3489566063dSJacob Faibussowitsch       PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2 * aircraft.nsteps, 2 * aircraft.nsteps, (void *)&aircraft, &aircraft.H));
3499566063dSJacob Faibussowitsch       PetscCall(MatShellSetOperation(aircraft.H, MATOP_MULT, (void (*)(void))MyMatMult));
3509566063dSJacob Faibussowitsch       PetscCall(MatSetOption(aircraft.H, MAT_SYMMETRIC, PETSC_TRUE));
3519566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao, aircraft.H, aircraft.H, MatrixFreeObjHessian, (void *)&aircraft));
352c4762a1bSJed Brown     } else {
3539566063dSJacob Faibussowitsch       PetscCall(MatCreateDense(MPI_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, 2 * aircraft.nsteps, 2 * aircraft.nsteps, NULL, &(aircraft.H)));
3549566063dSJacob Faibussowitsch       PetscCall(MatSetOption(aircraft.H, MAT_SYMMETRIC, PETSC_TRUE));
3559566063dSJacob Faibussowitsch       PetscCall(TaoSetHessian(tao, aircraft.H, aircraft.H, FormObjHessian, (void *)&aircraft));
356c4762a1bSJed Brown     }
357c4762a1bSJed Brown   }
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   /* Check for any TAO command line options */
3609566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao, &ksp));
361c4762a1bSJed Brown   if (ksp) {
3629566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
3639566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, PCNONE));
364c4762a1bSJed Brown   }
3659566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
366c4762a1bSJed Brown 
3679566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
3689566063dSJacob Faibussowitsch   PetscCall(VecView(P, PETSC_VIEWER_STDOUT_WORLD));
369c4762a1bSJed Brown 
370c4762a1bSJed Brown   /* Free TAO data structures */
3719566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
372c4762a1bSJed Brown 
373c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
374c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
375c4762a1bSJed Brown      are no longer needed.
376c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3779566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&aircraft.ts));
3789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&aircraft.A));
3799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aircraft.U));
3809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aircraft.V));
3819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aircraft.W));
3829566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&P));
3839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&PU));
3849566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&PL));
3859566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&aircraft.Jacp));
3869566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&aircraft.DRDU));
3879566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&aircraft.DRDP));
3889566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aircraft.Lambda[0]));
3899566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&aircraft.Mup[0]));
3909566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&P));
391c4762a1bSJed Brown   if (aircraft.eh) {
3929566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.Lambda2[0]));
3939566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.Mup2[0]));
3949566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.Dir));
3959566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.rhshp1[0]));
3969566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.rhshp2[0]));
3979566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.rhshp3[0]));
3989566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.rhshp4[0]));
3999566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.inthp1[0]));
4009566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.inthp2[0]));
4019566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.inthp3[0]));
4029566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&aircraft.inthp4[0]));
4039566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&aircraft.H));
404c4762a1bSJed Brown   }
4059566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
406b122ec5aSJacob Faibussowitsch   return 0;
407c4762a1bSJed Brown }
408c4762a1bSJed Brown 
409c4762a1bSJed Brown /*
410c4762a1bSJed Brown    FormObjFunctionGradient - Evaluates the function and corresponding gradient.
411c4762a1bSJed Brown 
412c4762a1bSJed Brown    Input Parameters:
413c4762a1bSJed Brown    tao - the Tao context
414c4762a1bSJed Brown    P   - the input vector
415a82e8c82SStefano Zampini    ctx - optional aircraft-defined context, as set by TaoSetObjectiveAndGradient()
416c4762a1bSJed Brown 
417c4762a1bSJed Brown    Output Parameters:
418c4762a1bSJed Brown    f   - the newly evaluated function
419c4762a1bSJed Brown    G   - the newly evaluated gradient
420c4762a1bSJed Brown */
421*9371c9d4SSatish Balay PetscErrorCode FormObjFunctionGradient(Tao tao, Vec P, PetscReal *f, Vec G, void *ctx) {
422c4762a1bSJed Brown   Aircraft           actx = (Aircraft)ctx;
423c4762a1bSJed Brown   TS                 ts   = actx->ts;
424c4762a1bSJed Brown   Vec                Q;
425c4762a1bSJed Brown   const PetscScalar *p, *q;
426c4762a1bSJed Brown   PetscScalar       *u, *v, *w;
427c4762a1bSJed Brown   PetscInt           i;
428c4762a1bSJed Brown 
429c4762a1bSJed Brown   PetscFunctionBeginUser;
4309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P, &p));
4319566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->V, &v));
4329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->W, &w));
433c4762a1bSJed Brown   for (i = 0; i < actx->nsteps; i++) {
434c4762a1bSJed Brown     v[i] = p[2 * i];
435c4762a1bSJed Brown     w[i] = p[2 * i + 1];
436c4762a1bSJed Brown   }
4379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P, &p));
4389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->V, &v));
4399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->W, &w));
440c4762a1bSJed Brown 
4419566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
4429566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
4439566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
4449566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, actx->ftime / actx->nsteps));
445c4762a1bSJed Brown 
446c4762a1bSJed Brown   /* reinitialize system state */
4479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->U, &u));
448c4762a1bSJed Brown   u[0] = 2.0;
449c4762a1bSJed Brown   u[1] = 0;
4509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->U, &u));
451c4762a1bSJed Brown 
452c4762a1bSJed Brown   /* reinitialize the integral value */
4539566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts, &Q));
4549566063dSJacob Faibussowitsch   PetscCall(VecSet(Q, 0.0));
455c4762a1bSJed Brown 
4569566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, actx->U));
457c4762a1bSJed Brown 
458c4762a1bSJed Brown   /* Reset initial conditions for the adjoint integration */
4599566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Lambda[0], 0.0));
4609566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Mup[0], 0.0));
4619566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, actx->Lambda, actx->Mup));
462c4762a1bSJed Brown 
4639566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
4649566063dSJacob Faibussowitsch   PetscCall(VecCopy(actx->Mup[0], G));
4659566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts, &Q));
4669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Q, &q));
467c4762a1bSJed Brown   *f = q[0];
4689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Q, &q));
469c4762a1bSJed Brown   PetscFunctionReturn(0);
470c4762a1bSJed Brown }
471c4762a1bSJed Brown 
472*9371c9d4SSatish Balay PetscErrorCode FormObjHessian(Tao tao, Vec P, Mat H, Mat Hpre, void *ctx) {
473c4762a1bSJed Brown   Aircraft           actx = (Aircraft)ctx;
474c4762a1bSJed Brown   const PetscScalar *p;
475c4762a1bSJed Brown   PetscScalar       *harr, *v, *w, one = 1.0;
476c4762a1bSJed Brown   PetscInt           ind[1];
477c4762a1bSJed Brown   PetscInt          *cols, i;
478c4762a1bSJed Brown   Vec                Dir;
479c4762a1bSJed Brown 
480c4762a1bSJed Brown   PetscFunctionBeginUser;
481c4762a1bSJed Brown   /* set up control parameters */
4829566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P, &p));
4839566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->V, &v));
4849566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->W, &w));
485c4762a1bSJed Brown   for (i = 0; i < actx->nsteps; i++) {
486c4762a1bSJed Brown     v[i] = p[2 * i];
487c4762a1bSJed Brown     w[i] = p[2 * i + 1];
488c4762a1bSJed Brown   }
4899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P, &p));
4909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->V, &v));
4919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->W, &w));
492c4762a1bSJed Brown 
4939566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * actx->nsteps, &harr));
4949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(2 * actx->nsteps, &cols));
495c4762a1bSJed Brown   for (i = 0; i < 2 * actx->nsteps; i++) cols[i] = i;
4969566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(P, &Dir));
497c4762a1bSJed Brown   for (i = 0; i < 2 * actx->nsteps; i++) {
498c4762a1bSJed Brown     ind[0] = i;
4999566063dSJacob Faibussowitsch     PetscCall(VecSet(Dir, 0.0));
5009566063dSJacob Faibussowitsch     PetscCall(VecSetValues(Dir, 1, ind, &one, INSERT_VALUES));
5019566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(Dir));
5029566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(Dir));
5039566063dSJacob Faibussowitsch     PetscCall(ComputeObjHessianWithSOA(Dir, harr, actx));
5049566063dSJacob Faibussowitsch     PetscCall(MatSetValues(H, 1, ind, 2 * actx->nsteps, cols, harr, INSERT_VALUES));
5059566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
5069566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
507c4762a1bSJed Brown     if (H != Hpre) {
5089566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(Hpre, MAT_FINAL_ASSEMBLY));
5099566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(Hpre, MAT_FINAL_ASSEMBLY));
510c4762a1bSJed Brown     }
511c4762a1bSJed Brown   }
5129566063dSJacob Faibussowitsch   PetscCall(PetscFree(cols));
5139566063dSJacob Faibussowitsch   PetscCall(PetscFree(harr));
5149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Dir));
515c4762a1bSJed Brown   PetscFunctionReturn(0);
516c4762a1bSJed Brown }
517c4762a1bSJed Brown 
518*9371c9d4SSatish Balay PetscErrorCode MatrixFreeObjHessian(Tao tao, Vec P, Mat H, Mat Hpre, void *ctx) {
519c4762a1bSJed Brown   Aircraft           actx = (Aircraft)ctx;
520c4762a1bSJed Brown   PetscScalar       *v, *w;
521c4762a1bSJed Brown   const PetscScalar *p;
522c4762a1bSJed Brown   PetscInt           i;
523c4762a1bSJed Brown 
524c4762a1bSJed Brown   PetscFunctionBegin;
5259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P, &p));
5269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->V, &v));
5279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->W, &w));
528c4762a1bSJed Brown   for (i = 0; i < actx->nsteps; i++) {
529c4762a1bSJed Brown     v[i] = p[2 * i];
530c4762a1bSJed Brown     w[i] = p[2 * i + 1];
531c4762a1bSJed Brown   }
5329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P, &p));
5339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->V, &v));
5349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->W, &w));
535c4762a1bSJed Brown   PetscFunctionReturn(0);
536c4762a1bSJed Brown }
537c4762a1bSJed Brown 
538*9371c9d4SSatish Balay PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y) {
539c4762a1bSJed Brown   PetscScalar *y;
540c4762a1bSJed Brown   void        *ptr;
541c4762a1bSJed Brown 
542c4762a1bSJed Brown   PetscFunctionBegin;
5439566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(H_shell, &ptr));
5449566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Y, &y));
5459566063dSJacob Faibussowitsch   PetscCall(ComputeObjHessianWithSOA(X, y, (Aircraft)ptr));
5469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Y, &y));
547c4762a1bSJed Brown   PetscFunctionReturn(0);
548c4762a1bSJed Brown }
549c4762a1bSJed Brown 
550*9371c9d4SSatish Balay PetscErrorCode ComputeObjHessianWithSOA(Vec Dir, PetscScalar arr[], Aircraft actx) {
551c4762a1bSJed Brown   TS                 ts = actx->ts;
552c4762a1bSJed Brown   const PetscScalar *z_ptr;
553c4762a1bSJed Brown   PetscScalar       *u;
554c4762a1bSJed Brown   Vec                Q;
555c4762a1bSJed Brown   PetscInt           i;
556c4762a1bSJed Brown 
557c4762a1bSJed Brown   PetscFunctionBeginUser;
558c4762a1bSJed Brown   /* Reset TSAdjoint so that AdjointSetUp will be called again */
5599566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
560c4762a1bSJed Brown 
5619566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
5629566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
5639566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
5649566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, actx->ftime / actx->nsteps));
5659566063dSJacob Faibussowitsch   PetscCall(TSSetCostHessianProducts(actx->ts, 1, actx->Lambda2, actx->Mup2, Dir));
566c4762a1bSJed Brown 
567c4762a1bSJed Brown   /* reinitialize system state */
5689566063dSJacob Faibussowitsch   PetscCall(VecGetArray(actx->U, &u));
569c4762a1bSJed Brown   u[0] = 2.0;
570c4762a1bSJed Brown   u[1] = 0;
5719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(actx->U, &u));
572c4762a1bSJed Brown 
573c4762a1bSJed Brown   /* reinitialize the integral value */
5749566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts, &Q));
5759566063dSJacob Faibussowitsch   PetscCall(VecSet(Q, 0.0));
576c4762a1bSJed Brown 
577c4762a1bSJed Brown   /* initialize tlm variable */
5789566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(actx->Jacp));
5799566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(actx->Jacp, MAT_FINAL_ASSEMBLY));
5809566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(actx->Jacp, MAT_FINAL_ASSEMBLY));
5819566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetForward(ts, actx->Jacp));
582c4762a1bSJed Brown 
5839566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, actx->U));
584c4762a1bSJed Brown 
585c4762a1bSJed Brown   /* Set terminal conditions for first- and second-order adjonts */
5869566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Lambda[0], 0.0));
5879566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Mup[0], 0.0));
5889566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Lambda2[0], 0.0));
5899566063dSJacob Faibussowitsch   PetscCall(VecSet(actx->Mup2[0], 0.0));
5909566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, actx->Lambda, actx->Mup));
591c4762a1bSJed Brown 
5929566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts, &Q));
593c4762a1bSJed Brown 
594c4762a1bSJed Brown   /* Reset initial conditions for the adjoint integration */
5959566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
596c4762a1bSJed Brown 
597c4762a1bSJed Brown   /* initial condition does not depend on p, so that lambda is not needed to assemble G */
5989566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(actx->Mup2[0], &z_ptr));
599c4762a1bSJed Brown   for (i = 0; i < 2 * actx->nsteps; i++) arr[i] = z_ptr[i];
6009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(actx->Mup2[0], &z_ptr));
601c4762a1bSJed Brown 
602c4762a1bSJed Brown   /* Disable second-order adjoint mode */
6039566063dSJacob Faibussowitsch   PetscCall(TSAdjointReset(ts));
6049566063dSJacob Faibussowitsch   PetscCall(TSAdjointResetForward(ts));
605c4762a1bSJed Brown   PetscFunctionReturn(0);
606c4762a1bSJed Brown }
607c4762a1bSJed Brown 
608c4762a1bSJed Brown /*TEST
609c4762a1bSJed Brown     build:
610c4762a1bSJed Brown       requires: !complex !single
611c4762a1bSJed Brown 
612c4762a1bSJed Brown     test:
613c4762a1bSJed Brown       args:  -ts_adapt_type none -ts_type rk -ts_rk_type 3 -viewer_binary_skip_info -tao_monitor -tao_gatol 1e-7
614c4762a1bSJed Brown 
615c4762a1bSJed Brown     test:
616c4762a1bSJed Brown       suffix: 2
617c4762a1bSJed 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
618c4762a1bSJed Brown 
619c4762a1bSJed Brown     test:
620c4762a1bSJed Brown       suffix: 3
621c4762a1bSJed 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
622c4762a1bSJed Brown TEST*/
623