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 24d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx) 25d71ae5a4SJacob Faibussowitsch { 26c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 27c4762a1bSJed Brown const PetscScalar *u, *v, *w; 28c4762a1bSJed Brown PetscScalar *f; 29c4762a1bSJed Brown PetscInt step; 30c4762a1bSJed Brown 31c4762a1bSJed Brown PetscFunctionBeginUser; 329566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &step)); 339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(actx->V, &v)); 359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(actx->W, &w)); 369566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 37c4762a1bSJed Brown f[0] = v[step] * PetscCosReal(w[step]); 38c4762a1bSJed Brown f[1] = v[step] * PetscSinReal(w[step]); 399566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(actx->V, &v)); 419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(actx->W, &w)); 429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec U, Mat A, void *ctx) 47d71ae5a4SJacob Faibussowitsch { 48c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 49c4762a1bSJed Brown const PetscScalar *u, *v, *w; 50c4762a1bSJed Brown PetscInt step, rows[2] = {0, 1}, rowcol[2]; 51c4762a1bSJed Brown PetscScalar Jp[2][2]; 52c4762a1bSJed Brown 53c4762a1bSJed Brown PetscFunctionBeginUser; 549566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(A)); 559566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &step)); 569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(actx->V, &v)); 589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(actx->W, &w)); 59c4762a1bSJed Brown 607b0e2f17SHong Zhang Jp[0][0] = PetscCosReal(w[step]); 617b0e2f17SHong Zhang Jp[0][1] = -v[step] * PetscSinReal(w[step]); 627b0e2f17SHong Zhang Jp[1][0] = PetscSinReal(w[step]); 637b0e2f17SHong Zhang Jp[1][1] = v[step] * PetscCosReal(w[step]); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 669566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(actx->V, &v)); 679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(actx->W, &w)); 68c4762a1bSJed Brown 697b0e2f17SHong Zhang rowcol[0] = 2 * step; 707b0e2f17SHong Zhang rowcol[1] = 2 * step + 1; 719566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, rows, 2, rowcol, &Jp[0][0], INSERT_VALUES)); 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 749566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 78d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 79d71ae5a4SJacob Faibussowitsch { 80c4762a1bSJed Brown PetscFunctionBeginUser; 813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 84d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSHessianProductUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 85d71ae5a4SJacob Faibussowitsch { 86c4762a1bSJed Brown PetscFunctionBeginUser; 873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSHessianProductPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 91d71ae5a4SJacob Faibussowitsch { 92c4762a1bSJed Brown PetscFunctionBeginUser; 933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSHessianProductPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 97d71ae5a4SJacob Faibussowitsch { 98c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 99c4762a1bSJed Brown const PetscScalar *v, *w, *vl, *vr, *u; 100c4762a1bSJed Brown PetscScalar *vhv; 101c4762a1bSJed Brown PetscScalar dJpdP[2][2][2] = {{{0}}}; 102c4762a1bSJed Brown PetscInt step, i, j, k; 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscFunctionBeginUser; 1059566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &step)); 1069566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1079566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(actx->V, &v)); 1089566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(actx->W, &w)); 1099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vl[0], &vl)); 1109566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vr, &vr)); 1119566063dSJacob Faibussowitsch PetscCall(VecSet(VHV[0], 0.0)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetArray(VHV[0], &vhv)); 113c4762a1bSJed Brown 1147b0e2f17SHong Zhang dJpdP[0][0][1] = -PetscSinReal(w[step]); 1157b0e2f17SHong Zhang dJpdP[0][1][0] = -PetscSinReal(w[step]); 1167b0e2f17SHong Zhang dJpdP[0][1][1] = -v[step] * PetscCosReal(w[step]); 1177b0e2f17SHong Zhang dJpdP[1][0][1] = PetscCosReal(w[step]); 1187b0e2f17SHong Zhang dJpdP[1][1][0] = PetscCosReal(w[step]); 1197b0e2f17SHong Zhang dJpdP[1][1][1] = -v[step] * PetscSinReal(w[step]); 120c4762a1bSJed Brown 121c4762a1bSJed Brown for (j = 0; j < 2; j++) { 1227b0e2f17SHong Zhang vhv[2 * step + j] = 0; 123c4762a1bSJed Brown for (k = 0; k < 2; k++) 1249371c9d4SSatish Balay for (i = 0; i < 2; i++) vhv[2 * step + j] += vl[i] * dJpdP[i][j][k] * vr[2 * step + k]; 125c4762a1bSJed Brown } 1269566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vl[0], &vl)); 1289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vr, &vr)); 1299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(VHV[0], &vhv)); 1303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* Vl in NULL,updates to VHV must be added */ 134d71ae5a4SJacob Faibussowitsch static PetscErrorCode IntegrandHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 135d71ae5a4SJacob Faibussowitsch { 136c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 137c4762a1bSJed Brown const PetscScalar *v, *w, *vr, *u; 138c4762a1bSJed Brown PetscScalar *vhv; 139c4762a1bSJed Brown PetscScalar dRudU[2][2] = {{0}}; 140c4762a1bSJed Brown PetscInt step, j, k; 141c4762a1bSJed Brown 142c4762a1bSJed Brown PetscFunctionBeginUser; 1439566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &step)); 1449566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1459566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(actx->V, &v)); 1469566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(actx->W, &w)); 1479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vr, &vr)); 1489566063dSJacob Faibussowitsch PetscCall(VecGetArray(VHV[0], &vhv)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown dRudU[0][0] = 2.0; 151c4762a1bSJed Brown dRudU[1][1] = 2.0; 152c4762a1bSJed Brown 153c4762a1bSJed Brown for (j = 0; j < 2; j++) { 154c4762a1bSJed Brown vhv[j] = 0; 1559371c9d4SSatish Balay for (k = 0; k < 2; k++) vhv[j] += dRudU[j][k] * vr[k]; 156c4762a1bSJed Brown } 1579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vr, &vr)); 1599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(VHV[0], &vhv)); 1603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163d71ae5a4SJacob Faibussowitsch static PetscErrorCode IntegrandHessianProductUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 164d71ae5a4SJacob Faibussowitsch { 165c4762a1bSJed Brown PetscFunctionBeginUser; 1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169d71ae5a4SJacob Faibussowitsch static PetscErrorCode IntegrandHessianProductPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 170d71ae5a4SJacob Faibussowitsch { 171c4762a1bSJed Brown PetscFunctionBeginUser; 1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173c4762a1bSJed Brown } 174c4762a1bSJed Brown 175d71ae5a4SJacob Faibussowitsch static PetscErrorCode IntegrandHessianProductPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 176d71ae5a4SJacob Faibussowitsch { 177c4762a1bSJed Brown PetscFunctionBeginUser; 1783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 179c4762a1bSJed Brown } 180c4762a1bSJed Brown 181d71ae5a4SJacob Faibussowitsch static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, void *ctx) 182d71ae5a4SJacob Faibussowitsch { 183c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 184c4762a1bSJed Brown PetscScalar *r; 185c4762a1bSJed Brown PetscReal dx, dy; 186c4762a1bSJed Brown const PetscScalar *u; 187c4762a1bSJed Brown 188c4762a1bSJed Brown PetscFunctionBegin; 1899566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1909566063dSJacob Faibussowitsch PetscCall(VecGetArray(R, &r)); 191c4762a1bSJed Brown dx = u[0] - actx->lv * t * PetscCosReal(actx->lw); 192c4762a1bSJed Brown dy = u[1] - actx->lv * t * PetscSinReal(actx->lw); 193c4762a1bSJed Brown r[0] = dx * dx + dy * dy; 1949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R, &r)); 1959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 197c4762a1bSJed Brown } 198c4762a1bSJed Brown 199d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, void *ctx) 200d71ae5a4SJacob Faibussowitsch { 201c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 202c4762a1bSJed Brown PetscScalar drdu[2][1]; 203c4762a1bSJed Brown const PetscScalar *u; 204c4762a1bSJed Brown PetscReal dx, dy; 205c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {0}; 206c4762a1bSJed Brown 207c4762a1bSJed Brown PetscFunctionBegin; 2089566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 209c4762a1bSJed Brown dx = u[0] - actx->lv * t * PetscCosReal(actx->lw); 210c4762a1bSJed Brown dy = u[1] - actx->lv * t * PetscSinReal(actx->lw); 211c4762a1bSJed Brown drdu[0][0] = 2. * dx; 212c4762a1bSJed Brown drdu[1][0] = 2. * dy; 2139566063dSJacob Faibussowitsch PetscCall(MatSetValues(DRDU, 2, row, 1, col, &drdu[0][0], INSERT_VALUES)); 2149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 2159566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY)); 2169566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY)); 2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 220d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, void *ctx) 221d71ae5a4SJacob Faibussowitsch { 222c4762a1bSJed Brown PetscFunctionBegin; 2239566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(DRDP)); 2243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 225c4762a1bSJed Brown } 226c4762a1bSJed Brown 227d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 228d71ae5a4SJacob Faibussowitsch { 229c4762a1bSJed Brown Vec P, PL, PU; 230c4762a1bSJed Brown struct _n_aircraft aircraft; 231c4762a1bSJed Brown PetscMPIInt size; 232c4762a1bSJed Brown Tao tao; 233c4762a1bSJed Brown KSP ksp; 234c4762a1bSJed Brown PC pc; 235c4762a1bSJed Brown PetscScalar *u, *p; 236c4762a1bSJed Brown PetscInt i; 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* Initialize program */ 239327415f7SBarry Smith PetscFunctionBeginUser; 2409566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, NULL)); 2419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2423c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 243c4762a1bSJed Brown 244c4762a1bSJed Brown /* Parameter settings */ 245c4762a1bSJed Brown aircraft.ftime = 1.; /* time interval in hour */ 246c4762a1bSJed Brown aircraft.nsteps = 10; /* number of steps */ 247c4762a1bSJed Brown aircraft.lv = 2.0; /* leader speed in kmph */ 248c4762a1bSJed Brown aircraft.lw = PETSC_PI / 4.; /* leader heading angle */ 249c4762a1bSJed Brown 2509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-ftime", &aircraft.ftime, NULL)); 2519566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nsteps", &aircraft.nsteps, NULL)); 2529566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-matrixfree", &aircraft.mf)); 2539566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-exacthessian", &aircraft.eh)); 254c4762a1bSJed Brown 255c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 2569566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 2579566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBQNLS)); 258c4762a1bSJed Brown 259c4762a1bSJed Brown /* Create necessary matrix and vectors, solve same ODE on every process */ 2609566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &aircraft.A)); 2619566063dSJacob Faibussowitsch PetscCall(MatSetSizes(aircraft.A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 2629566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(aircraft.A)); 2639566063dSJacob Faibussowitsch PetscCall(MatSetUp(aircraft.A)); 26426cec326SBarry Smith /* this is to set explicit zeros along the diagonal of the matrix */ 2659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(aircraft.A, MAT_FINAL_ASSEMBLY)); 2669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(aircraft.A, MAT_FINAL_ASSEMBLY)); 2679566063dSJacob Faibussowitsch PetscCall(MatShift(aircraft.A, 1)); 2689566063dSJacob Faibussowitsch PetscCall(MatShift(aircraft.A, -1)); 269c4762a1bSJed Brown 2709566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &aircraft.Jacp)); 2719566063dSJacob Faibussowitsch PetscCall(MatSetSizes(aircraft.Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 2 * aircraft.nsteps)); 2729566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(aircraft.Jacp)); 2739566063dSJacob Faibussowitsch PetscCall(MatSetUp(aircraft.Jacp)); 27426cec326SBarry Smith PetscCall(MatSetOption(aircraft.Jacp, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 2759566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2 * aircraft.nsteps, 1, NULL, &aircraft.DRDP)); 2769566063dSJacob Faibussowitsch PetscCall(MatSetUp(aircraft.DRDP)); 2779566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, &aircraft.DRDU)); 2789566063dSJacob Faibussowitsch PetscCall(MatSetUp(aircraft.DRDU)); 279c4762a1bSJed Brown 280c4762a1bSJed Brown /* Create timestepping solver context */ 2819566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &aircraft.ts)); 2829566063dSJacob Faibussowitsch PetscCall(TSSetType(aircraft.ts, TSRK)); 2839566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(aircraft.ts, NULL, RHSFunction, &aircraft)); 2849566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(aircraft.ts, aircraft.A, aircraft.A, TSComputeRHSJacobianConstant, &aircraft)); 2859566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(aircraft.ts, aircraft.Jacp, RHSJacobianP, &aircraft)); 2869566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(aircraft.ts, TS_EXACTFINALTIME_MATCHSTEP)); 2879566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(aircraft.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 288c4762a1bSJed Brown 289c4762a1bSJed Brown /* Set initial conditions */ 2909566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.A, &aircraft.U, NULL)); 2919566063dSJacob Faibussowitsch PetscCall(TSSetSolution(aircraft.ts, aircraft.U)); 2929566063dSJacob Faibussowitsch PetscCall(VecGetArray(aircraft.U, &u)); 293c4762a1bSJed Brown u[0] = 1.5; 294c4762a1bSJed Brown u[1] = 0; 2959566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(aircraft.U, &u)); 2969566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &aircraft.V)); 2979566063dSJacob Faibussowitsch PetscCall(VecSetSizes(aircraft.V, PETSC_DECIDE, aircraft.nsteps)); 2989566063dSJacob Faibussowitsch PetscCall(VecSetUp(aircraft.V)); 2999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(aircraft.V, &aircraft.W)); 3009566063dSJacob Faibussowitsch PetscCall(VecSet(aircraft.V, 1.)); 3019566063dSJacob Faibussowitsch PetscCall(VecSet(aircraft.W, PETSC_PI / 4.)); 302c4762a1bSJed Brown 303c4762a1bSJed Brown /* Save trajectory of solution so that TSAdjointSolve() may be used */ 3049566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(aircraft.ts)); 305c4762a1bSJed Brown 306c4762a1bSJed Brown /* Set sensitivity context */ 3079566063dSJacob Faibussowitsch PetscCall(TSCreateQuadratureTS(aircraft.ts, PETSC_FALSE, &aircraft.quadts)); 3088434afd1SBarry Smith PetscCall(TSSetRHSFunction(aircraft.quadts, NULL, (TSRHSFunctionFn *)CostIntegrand, &aircraft)); 3099566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(aircraft.quadts, aircraft.DRDU, aircraft.DRDU, (TSRHSJacobian)DRDUJacobianTranspose, &aircraft)); 3108434afd1SBarry Smith PetscCall(TSSetRHSJacobianP(aircraft.quadts, aircraft.DRDP, (TSRHSJacobianPFn *)DRDPJacobianTranspose, &aircraft)); 3119566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.A, &aircraft.Lambda[0], NULL)); 3129566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.Jacp, &aircraft.Mup[0], NULL)); 313c4762a1bSJed Brown if (aircraft.eh) { 3149566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.A, &aircraft.rhshp1[0], NULL)); 3159566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.A, &aircraft.rhshp2[0], NULL)); 3169566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.Jacp, &aircraft.rhshp3[0], NULL)); 3179566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.Jacp, &aircraft.rhshp4[0], NULL)); 3189566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.DRDU, &aircraft.inthp1[0], NULL)); 3199566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.DRDU, &aircraft.inthp2[0], NULL)); 3209566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.DRDP, &aircraft.inthp3[0], NULL)); 3219566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.DRDP, &aircraft.inthp4[0], NULL)); 3229566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.Jacp, &aircraft.Dir, NULL)); 3239566063dSJacob Faibussowitsch PetscCall(TSSetRHSHessianProduct(aircraft.ts, aircraft.rhshp1, RHSHessianProductUU, aircraft.rhshp2, RHSHessianProductUP, aircraft.rhshp3, RHSHessianProductPU, aircraft.rhshp4, RHSHessianProductPP, &aircraft)); 3249566063dSJacob Faibussowitsch PetscCall(TSSetRHSHessianProduct(aircraft.quadts, aircraft.inthp1, IntegrandHessianProductUU, aircraft.inthp2, IntegrandHessianProductUP, aircraft.inthp3, IntegrandHessianProductPU, aircraft.inthp4, IntegrandHessianProductPP, &aircraft)); 3259566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.A, &aircraft.Lambda2[0], NULL)); 3269566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.Jacp, &aircraft.Mup2[0], NULL)); 327c4762a1bSJed Brown } 3289566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(aircraft.ts)); 3299566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(aircraft.ts, aircraft.ftime)); 3309566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(aircraft.ts, aircraft.ftime / aircraft.nsteps)); 331c4762a1bSJed Brown 332c4762a1bSJed Brown /* Set initial solution guess */ 3339566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(aircraft.Jacp, &P, NULL)); 3349566063dSJacob Faibussowitsch PetscCall(VecGetArray(P, &p)); 335c4762a1bSJed Brown for (i = 0; i < aircraft.nsteps; i++) { 336c4762a1bSJed Brown p[2 * i] = 2.0; 337c4762a1bSJed Brown p[2 * i + 1] = PETSC_PI / 2.0; 338c4762a1bSJed Brown } 3399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(P, &p)); 3409566063dSJacob Faibussowitsch PetscCall(VecDuplicate(P, &PU)); 3419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(P, &PL)); 3429566063dSJacob Faibussowitsch PetscCall(VecGetArray(PU, &p)); 343c4762a1bSJed Brown for (i = 0; i < aircraft.nsteps; i++) { 344c4762a1bSJed Brown p[2 * i] = 2.0; 345c4762a1bSJed Brown p[2 * i + 1] = PETSC_PI; 346c4762a1bSJed Brown } 3479566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(PU, &p)); 3489566063dSJacob Faibussowitsch PetscCall(VecGetArray(PL, &p)); 349c4762a1bSJed Brown for (i = 0; i < aircraft.nsteps; i++) { 350c4762a1bSJed Brown p[2 * i] = 0.0; 351c4762a1bSJed Brown p[2 * i + 1] = -PETSC_PI; 352c4762a1bSJed Brown } 3539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(PL, &p)); 354c4762a1bSJed Brown 3559566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, P)); 3569566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao, PL, PU)); 357c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 3589566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormObjFunctionGradient, (void *)&aircraft)); 359c4762a1bSJed Brown 360c4762a1bSJed Brown if (aircraft.eh) { 361c4762a1bSJed Brown if (aircraft.mf) { 3629566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2 * aircraft.nsteps, 2 * aircraft.nsteps, (void *)&aircraft, &aircraft.H)); 3639566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(aircraft.H, MATOP_MULT, (void (*)(void))MyMatMult)); 3649566063dSJacob Faibussowitsch PetscCall(MatSetOption(aircraft.H, MAT_SYMMETRIC, PETSC_TRUE)); 3659566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, aircraft.H, aircraft.H, MatrixFreeObjHessian, (void *)&aircraft)); 366c4762a1bSJed Brown } else { 367*f4f49eeaSPierre Jolivet PetscCall(MatCreateDense(MPI_COMM_WORLD, PETSC_DETERMINE, PETSC_DETERMINE, 2 * aircraft.nsteps, 2 * aircraft.nsteps, NULL, &aircraft.H)); 3689566063dSJacob Faibussowitsch PetscCall(MatSetOption(aircraft.H, MAT_SYMMETRIC, PETSC_TRUE)); 3699566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, aircraft.H, aircraft.H, FormObjHessian, (void *)&aircraft)); 370c4762a1bSJed Brown } 371c4762a1bSJed Brown } 372c4762a1bSJed Brown 373c4762a1bSJed Brown /* Check for any TAO command line options */ 3749566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 375c4762a1bSJed Brown if (ksp) { 3769566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 3779566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 378c4762a1bSJed Brown } 3799566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 380c4762a1bSJed Brown 3819566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 3829566063dSJacob Faibussowitsch PetscCall(VecView(P, PETSC_VIEWER_STDOUT_WORLD)); 383c4762a1bSJed Brown 384c4762a1bSJed Brown /* Free TAO data structures */ 3859566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 386c4762a1bSJed Brown 387c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 388c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 389c4762a1bSJed Brown are no longer needed. 390c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3919566063dSJacob Faibussowitsch PetscCall(TSDestroy(&aircraft.ts)); 3929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aircraft.A)); 3939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.U)); 3949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.V)); 3959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.W)); 3969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&P)); 3979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&PU)); 3989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&PL)); 3999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aircraft.Jacp)); 4009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aircraft.DRDU)); 4019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aircraft.DRDP)); 4029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.Lambda[0])); 4039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.Mup[0])); 4049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&P)); 405c4762a1bSJed Brown if (aircraft.eh) { 4069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.Lambda2[0])); 4079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.Mup2[0])); 4089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.Dir)); 4099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.rhshp1[0])); 4109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.rhshp2[0])); 4119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.rhshp3[0])); 4129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.rhshp4[0])); 4139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.inthp1[0])); 4149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.inthp2[0])); 4159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.inthp3[0])); 4169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&aircraft.inthp4[0])); 4179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&aircraft.H)); 418c4762a1bSJed Brown } 4199566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 420b122ec5aSJacob Faibussowitsch return 0; 421c4762a1bSJed Brown } 422c4762a1bSJed Brown 423c4762a1bSJed Brown /* 424c4762a1bSJed Brown FormObjFunctionGradient - Evaluates the function and corresponding gradient. 425c4762a1bSJed Brown 426c4762a1bSJed Brown Input Parameters: 427c4762a1bSJed Brown tao - the Tao context 428c4762a1bSJed Brown P - the input vector 429a82e8c82SStefano Zampini ctx - optional aircraft-defined context, as set by TaoSetObjectiveAndGradient() 430c4762a1bSJed Brown 431c4762a1bSJed Brown Output Parameters: 432c4762a1bSJed Brown f - the newly evaluated function 433c4762a1bSJed Brown G - the newly evaluated gradient 434c4762a1bSJed Brown */ 435d71ae5a4SJacob Faibussowitsch PetscErrorCode FormObjFunctionGradient(Tao tao, Vec P, PetscReal *f, Vec G, void *ctx) 436d71ae5a4SJacob Faibussowitsch { 437c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 438c4762a1bSJed Brown TS ts = actx->ts; 439c4762a1bSJed Brown Vec Q; 440c4762a1bSJed Brown const PetscScalar *p, *q; 441c4762a1bSJed Brown PetscScalar *u, *v, *w; 442c4762a1bSJed Brown PetscInt i; 443c4762a1bSJed Brown 444c4762a1bSJed Brown PetscFunctionBeginUser; 4459566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, &p)); 4469566063dSJacob Faibussowitsch PetscCall(VecGetArray(actx->V, &v)); 4479566063dSJacob Faibussowitsch PetscCall(VecGetArray(actx->W, &w)); 448c4762a1bSJed Brown for (i = 0; i < actx->nsteps; i++) { 449c4762a1bSJed Brown v[i] = p[2 * i]; 450c4762a1bSJed Brown w[i] = p[2 * i + 1]; 451c4762a1bSJed Brown } 4529566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, &p)); 4539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(actx->V, &v)); 4549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(actx->W, &w)); 455c4762a1bSJed Brown 4569566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 4579566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 4589566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 4599566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, actx->ftime / actx->nsteps)); 460c4762a1bSJed Brown 461c4762a1bSJed Brown /* reinitialize system state */ 4629566063dSJacob Faibussowitsch PetscCall(VecGetArray(actx->U, &u)); 463c4762a1bSJed Brown u[0] = 2.0; 464c4762a1bSJed Brown u[1] = 0; 4659566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(actx->U, &u)); 466c4762a1bSJed Brown 467c4762a1bSJed Brown /* reinitialize the integral value */ 4689566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &Q)); 4699566063dSJacob Faibussowitsch PetscCall(VecSet(Q, 0.0)); 470c4762a1bSJed Brown 4719566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, actx->U)); 472c4762a1bSJed Brown 473c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 4749566063dSJacob Faibussowitsch PetscCall(VecSet(actx->Lambda[0], 0.0)); 4759566063dSJacob Faibussowitsch PetscCall(VecSet(actx->Mup[0], 0.0)); 4769566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, actx->Lambda, actx->Mup)); 477c4762a1bSJed Brown 4789566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 4799566063dSJacob Faibussowitsch PetscCall(VecCopy(actx->Mup[0], G)); 4809566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &Q)); 4819566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Q, &q)); 482c4762a1bSJed Brown *f = q[0]; 4839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Q, &q)); 4843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 485c4762a1bSJed Brown } 486c4762a1bSJed Brown 487d71ae5a4SJacob Faibussowitsch PetscErrorCode FormObjHessian(Tao tao, Vec P, Mat H, Mat Hpre, void *ctx) 488d71ae5a4SJacob Faibussowitsch { 489c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 490c4762a1bSJed Brown const PetscScalar *p; 491c4762a1bSJed Brown PetscScalar *harr, *v, *w, one = 1.0; 492c4762a1bSJed Brown PetscInt ind[1]; 493c4762a1bSJed Brown PetscInt *cols, i; 494c4762a1bSJed Brown Vec Dir; 495c4762a1bSJed Brown 496c4762a1bSJed Brown PetscFunctionBeginUser; 497c4762a1bSJed Brown /* set up control parameters */ 4989566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, &p)); 4999566063dSJacob Faibussowitsch PetscCall(VecGetArray(actx->V, &v)); 5009566063dSJacob Faibussowitsch PetscCall(VecGetArray(actx->W, &w)); 501c4762a1bSJed Brown for (i = 0; i < actx->nsteps; i++) { 502c4762a1bSJed Brown v[i] = p[2 * i]; 503c4762a1bSJed Brown w[i] = p[2 * i + 1]; 504c4762a1bSJed Brown } 5059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, &p)); 5069566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(actx->V, &v)); 5079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(actx->W, &w)); 508c4762a1bSJed Brown 5099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * actx->nsteps, &harr)); 5109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * actx->nsteps, &cols)); 511c4762a1bSJed Brown for (i = 0; i < 2 * actx->nsteps; i++) cols[i] = i; 5129566063dSJacob Faibussowitsch PetscCall(VecDuplicate(P, &Dir)); 513c4762a1bSJed Brown for (i = 0; i < 2 * actx->nsteps; i++) { 514c4762a1bSJed Brown ind[0] = i; 5159566063dSJacob Faibussowitsch PetscCall(VecSet(Dir, 0.0)); 5169566063dSJacob Faibussowitsch PetscCall(VecSetValues(Dir, 1, ind, &one, INSERT_VALUES)); 5179566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Dir)); 5189566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Dir)); 5199566063dSJacob Faibussowitsch PetscCall(ComputeObjHessianWithSOA(Dir, harr, actx)); 5209566063dSJacob Faibussowitsch PetscCall(MatSetValues(H, 1, ind, 2 * actx->nsteps, cols, harr, INSERT_VALUES)); 5219566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY)); 5229566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY)); 523c4762a1bSJed Brown if (H != Hpre) { 5249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Hpre, MAT_FINAL_ASSEMBLY)); 5259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Hpre, MAT_FINAL_ASSEMBLY)); 526c4762a1bSJed Brown } 527c4762a1bSJed Brown } 5289566063dSJacob Faibussowitsch PetscCall(PetscFree(cols)); 5299566063dSJacob Faibussowitsch PetscCall(PetscFree(harr)); 5309566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Dir)); 5313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 532c4762a1bSJed Brown } 533c4762a1bSJed Brown 534d71ae5a4SJacob Faibussowitsch PetscErrorCode MatrixFreeObjHessian(Tao tao, Vec P, Mat H, Mat Hpre, void *ctx) 535d71ae5a4SJacob Faibussowitsch { 536c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 537c4762a1bSJed Brown PetscScalar *v, *w; 538c4762a1bSJed Brown const PetscScalar *p; 539c4762a1bSJed Brown PetscInt i; 540c4762a1bSJed Brown 541c4762a1bSJed Brown PetscFunctionBegin; 5429566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, &p)); 5439566063dSJacob Faibussowitsch PetscCall(VecGetArray(actx->V, &v)); 5449566063dSJacob Faibussowitsch PetscCall(VecGetArray(actx->W, &w)); 545c4762a1bSJed Brown for (i = 0; i < actx->nsteps; i++) { 546c4762a1bSJed Brown v[i] = p[2 * i]; 547c4762a1bSJed Brown w[i] = p[2 * i + 1]; 548c4762a1bSJed Brown } 5499566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, &p)); 5509566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(actx->V, &v)); 5519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(actx->W, &w)); 5523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 553c4762a1bSJed Brown } 554c4762a1bSJed Brown 555d71ae5a4SJacob Faibussowitsch PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y) 556d71ae5a4SJacob Faibussowitsch { 557c4762a1bSJed Brown PetscScalar *y; 558c4762a1bSJed Brown void *ptr; 559c4762a1bSJed Brown 560c4762a1bSJed Brown PetscFunctionBegin; 5619566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(H_shell, &ptr)); 5629566063dSJacob Faibussowitsch PetscCall(VecGetArray(Y, &y)); 5639566063dSJacob Faibussowitsch PetscCall(ComputeObjHessianWithSOA(X, y, (Aircraft)ptr)); 5649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Y, &y)); 5653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 566c4762a1bSJed Brown } 567c4762a1bSJed Brown 568d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeObjHessianWithSOA(Vec Dir, PetscScalar arr[], Aircraft actx) 569d71ae5a4SJacob Faibussowitsch { 570c4762a1bSJed Brown TS ts = actx->ts; 571c4762a1bSJed Brown const PetscScalar *z_ptr; 572c4762a1bSJed Brown PetscScalar *u; 573c4762a1bSJed Brown Vec Q; 574c4762a1bSJed Brown PetscInt i; 575c4762a1bSJed Brown 576c4762a1bSJed Brown PetscFunctionBeginUser; 577c4762a1bSJed Brown /* Reset TSAdjoint so that AdjointSetUp will be called again */ 5789566063dSJacob Faibussowitsch PetscCall(TSAdjointReset(ts)); 579c4762a1bSJed Brown 5809566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 5819566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 5829566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 5839566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, actx->ftime / actx->nsteps)); 5849566063dSJacob Faibussowitsch PetscCall(TSSetCostHessianProducts(actx->ts, 1, actx->Lambda2, actx->Mup2, Dir)); 585c4762a1bSJed Brown 586c4762a1bSJed Brown /* reinitialize system state */ 5879566063dSJacob Faibussowitsch PetscCall(VecGetArray(actx->U, &u)); 588c4762a1bSJed Brown u[0] = 2.0; 589c4762a1bSJed Brown u[1] = 0; 5909566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(actx->U, &u)); 591c4762a1bSJed Brown 592c4762a1bSJed Brown /* reinitialize the integral value */ 5939566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &Q)); 5949566063dSJacob Faibussowitsch PetscCall(VecSet(Q, 0.0)); 595c4762a1bSJed Brown 596c4762a1bSJed Brown /* initialize tlm variable */ 5979566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(actx->Jacp)); 5989566063dSJacob Faibussowitsch PetscCall(TSAdjointSetForward(ts, actx->Jacp)); 599c4762a1bSJed Brown 6009566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, actx->U)); 601c4762a1bSJed Brown 602c4762a1bSJed Brown /* Set terminal conditions for first- and second-order adjonts */ 6039566063dSJacob Faibussowitsch PetscCall(VecSet(actx->Lambda[0], 0.0)); 6049566063dSJacob Faibussowitsch PetscCall(VecSet(actx->Mup[0], 0.0)); 6059566063dSJacob Faibussowitsch PetscCall(VecSet(actx->Lambda2[0], 0.0)); 6069566063dSJacob Faibussowitsch PetscCall(VecSet(actx->Mup2[0], 0.0)); 6079566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, actx->Lambda, actx->Mup)); 608c4762a1bSJed Brown 6099566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &Q)); 610c4762a1bSJed Brown 611c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 6129566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 613c4762a1bSJed Brown 614c4762a1bSJed Brown /* initial condition does not depend on p, so that lambda is not needed to assemble G */ 6159566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(actx->Mup2[0], &z_ptr)); 616c4762a1bSJed Brown for (i = 0; i < 2 * actx->nsteps; i++) arr[i] = z_ptr[i]; 6179566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(actx->Mup2[0], &z_ptr)); 618c4762a1bSJed Brown 619c4762a1bSJed Brown /* Disable second-order adjoint mode */ 6209566063dSJacob Faibussowitsch PetscCall(TSAdjointReset(ts)); 6219566063dSJacob Faibussowitsch PetscCall(TSAdjointResetForward(ts)); 6223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 623c4762a1bSJed Brown } 624c4762a1bSJed Brown 625c4762a1bSJed Brown /*TEST 626c4762a1bSJed Brown build: 627c4762a1bSJed Brown requires: !complex !single 628c4762a1bSJed Brown 629c4762a1bSJed Brown test: 630c4762a1bSJed Brown args: -ts_adapt_type none -ts_type rk -ts_rk_type 3 -viewer_binary_skip_info -tao_monitor -tao_gatol 1e-7 631c4762a1bSJed Brown 632c4762a1bSJed Brown test: 633c4762a1bSJed Brown suffix: 2 634c4762a1bSJed 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 635c4762a1bSJed Brown 636c4762a1bSJed Brown test: 637c4762a1bSJed Brown suffix: 3 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 -matrixfree 639c4762a1bSJed Brown TEST*/ 640