1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation.\n\ 3c4762a1bSJed Brown Input parameters include:\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown /* ------------------------------------------------------------------------ 6c4762a1bSJed Brown 7c4762a1bSJed Brown Notes: 8c4762a1bSJed Brown This code demonstrates how to solve a DAE-constrained optimization problem with TAO, TSAdjoint and TS. 9c4762a1bSJed Brown The nonlinear problem is written in a DAE equivalent form. 10c4762a1bSJed Brown The objective is to minimize the difference between observation and model prediction by finding an optimal value for parameter \mu. 11c4762a1bSJed Brown The gradient is computed with the discrete adjoint of an implicit theta method, see ex20adj.c for details. 12c4762a1bSJed Brown ------------------------------------------------------------------------- */ 13c4762a1bSJed Brown #include <petsctao.h> 14c4762a1bSJed Brown #include <petscts.h> 15c4762a1bSJed Brown 16c4762a1bSJed Brown typedef struct _n_User *User; 17c4762a1bSJed Brown struct _n_User { 18c4762a1bSJed Brown TS ts; 19c4762a1bSJed Brown PetscReal mu; 20c4762a1bSJed Brown PetscReal next_output; 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* Sensitivity analysis support */ 23c4762a1bSJed Brown PetscReal ftime; 24c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 25c4762a1bSJed Brown Mat Jacp; /* JacobianP matrix */ 26c4762a1bSJed Brown Mat H; /* Hessian matrix for optimization */ 27c4762a1bSJed Brown Vec U, Lambda[1], Mup[1]; /* adjoint variables */ 28c4762a1bSJed Brown Vec Lambda2[1], Mup2[1]; /* second-order adjoint variables */ 29c4762a1bSJed Brown Vec Ihp1[1]; /* working space for Hessian evaluations */ 30c4762a1bSJed Brown Vec Ihp2[1]; /* working space for Hessian evaluations */ 31c4762a1bSJed Brown Vec Ihp3[1]; /* working space for Hessian evaluations */ 32c4762a1bSJed Brown Vec Ihp4[1]; /* working space for Hessian evaluations */ 33c4762a1bSJed Brown Vec Dir; /* direction vector */ 34c4762a1bSJed Brown PetscReal ob[2]; /* observation used by the cost function */ 35c4762a1bSJed Brown PetscBool implicitform; /* implicit ODE? */ 36c4762a1bSJed Brown }; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 39c4762a1bSJed Brown PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *); 40c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec, PetscScalar[], User); 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE -------------------- */ 43c4762a1bSJed Brown 44d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx) 45d71ae5a4SJacob Faibussowitsch { 46c4762a1bSJed Brown User user = (User)ctx; 47c4762a1bSJed Brown PetscScalar *f; 48c4762a1bSJed Brown const PetscScalar *u; 49c4762a1bSJed Brown 50c4762a1bSJed Brown PetscFunctionBeginUser; 519566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 529566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 53c4762a1bSJed Brown f[0] = u[1]; 54c4762a1bSJed Brown f[1] = user->mu * ((1. - u[0] * u[0]) * u[1] - u[0]); 559566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 57*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx) 61d71ae5a4SJacob Faibussowitsch { 62c4762a1bSJed Brown User user = (User)ctx; 63c4762a1bSJed Brown PetscReal mu = user->mu; 64c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 65c4762a1bSJed Brown PetscScalar J[2][2]; 66c4762a1bSJed Brown const PetscScalar *u; 67c4762a1bSJed Brown 68c4762a1bSJed Brown PetscFunctionBeginUser; 699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 70c4762a1bSJed Brown 71c4762a1bSJed Brown J[0][0] = 0; 72c4762a1bSJed Brown J[1][0] = -mu * (2.0 * u[1] * u[0] + 1.); 73c4762a1bSJed Brown J[0][1] = 1.0; 74c4762a1bSJed Brown J[1][1] = mu * (1.0 - u[0] * u[0]); 759566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 769566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 779566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 78c4762a1bSJed Brown if (B && A != B) { 799566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 809566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 84*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 85c4762a1bSJed Brown } 86c4762a1bSJed Brown 87d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 88d71ae5a4SJacob Faibussowitsch { 89c4762a1bSJed Brown const PetscScalar *vl, *vr, *u; 90c4762a1bSJed Brown PetscScalar *vhv; 91c4762a1bSJed Brown PetscScalar dJdU[2][2][2] = {{{0}}}; 92c4762a1bSJed Brown PetscInt i, j, k; 93c4762a1bSJed Brown User user = (User)ctx; 94c4762a1bSJed Brown 95c4762a1bSJed Brown PetscFunctionBeginUser; 969566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 979566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vl[0], &vl)); 989566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vr, &vr)); 999566063dSJacob Faibussowitsch PetscCall(VecGetArray(VHV[0], &vhv)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown dJdU[1][0][0] = -2. * user->mu * u[1]; 102c4762a1bSJed Brown dJdU[1][1][0] = -2. * user->mu * u[0]; 103c4762a1bSJed Brown dJdU[1][0][1] = -2. * user->mu * u[0]; 104c4762a1bSJed Brown for (j = 0; j < 2; j++) { 105c4762a1bSJed Brown vhv[j] = 0; 106c4762a1bSJed Brown for (k = 0; k < 2; k++) 1079371c9d4SSatish Balay for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k]; 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 1109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vl[0], &vl)); 1129566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vr, &vr)); 1139566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(VHV[0], &vhv)); 114*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown 117d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSHessianProductUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 118d71ae5a4SJacob Faibussowitsch { 119c4762a1bSJed Brown const PetscScalar *vl, *vr, *u; 120c4762a1bSJed Brown PetscScalar *vhv; 121c4762a1bSJed Brown PetscScalar dJdP[2][2][1] = {{{0}}}; 122c4762a1bSJed Brown PetscInt i, j, k; 123c4762a1bSJed Brown 124c4762a1bSJed Brown PetscFunctionBeginUser; 1259566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1269566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vl[0], &vl)); 1279566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vr, &vr)); 1289566063dSJacob Faibussowitsch PetscCall(VecGetArray(VHV[0], &vhv)); 129c4762a1bSJed Brown 130c4762a1bSJed Brown dJdP[1][0][0] = -(1. + 2. * u[0] * u[1]); 131c4762a1bSJed Brown dJdP[1][1][0] = 1. - u[0] * u[0]; 132c4762a1bSJed Brown for (j = 0; j < 2; j++) { 133c4762a1bSJed Brown vhv[j] = 0; 134c4762a1bSJed Brown for (k = 0; k < 1; k++) 1359371c9d4SSatish Balay for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdP[i][j][k] * vr[k]; 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 1389566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1399566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vl[0], &vl)); 1409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vr, &vr)); 1419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(VHV[0], &vhv)); 142*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSHessianProductPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 146d71ae5a4SJacob Faibussowitsch { 147c4762a1bSJed Brown const PetscScalar *vl, *vr, *u; 148c4762a1bSJed Brown PetscScalar *vhv; 149c4762a1bSJed Brown PetscScalar dJdU[2][1][2] = {{{0}}}; 150c4762a1bSJed Brown PetscInt i, j, k; 151c4762a1bSJed Brown 152c4762a1bSJed Brown PetscFunctionBeginUser; 1539566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1549566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vl[0], &vl)); 1559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vr, &vr)); 1569566063dSJacob Faibussowitsch PetscCall(VecGetArray(VHV[0], &vhv)); 157c4762a1bSJed Brown 158c4762a1bSJed Brown dJdU[1][0][0] = -1. - 2. * u[1] * u[0]; 159c4762a1bSJed Brown dJdU[1][0][1] = 1. - u[0] * u[0]; 160c4762a1bSJed Brown for (j = 0; j < 1; j++) { 161c4762a1bSJed Brown vhv[j] = 0; 162c4762a1bSJed Brown for (k = 0; k < 2; k++) 1639371c9d4SSatish Balay for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k]; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 1669566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vl[0], &vl)); 1689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vr, &vr)); 1699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(VHV[0], &vhv)); 170*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSHessianProductPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 174d71ae5a4SJacob Faibussowitsch { 175c4762a1bSJed Brown PetscFunctionBeginUser; 176*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE -------------------- */ 180c4762a1bSJed Brown 181d71ae5a4SJacob Faibussowitsch static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) 182d71ae5a4SJacob Faibussowitsch { 183c4762a1bSJed Brown User user = (User)ctx; 184c4762a1bSJed Brown PetscScalar *f; 185c4762a1bSJed Brown const PetscScalar *u, *udot; 186c4762a1bSJed Brown 187c4762a1bSJed Brown PetscFunctionBeginUser; 1889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1899566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 1909566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 191c4762a1bSJed Brown 192c4762a1bSJed Brown f[0] = udot[0] - u[1]; 193c4762a1bSJed Brown f[1] = udot[1] - user->mu * ((1.0 - u[0] * u[0]) * u[1] - u[0]); 194c4762a1bSJed Brown 1959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 1979566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 198*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201d71ae5a4SJacob Faibussowitsch static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) 202d71ae5a4SJacob Faibussowitsch { 203c4762a1bSJed Brown User user = (User)ctx; 204c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 205c4762a1bSJed Brown PetscScalar J[2][2]; 206c4762a1bSJed Brown const PetscScalar *u; 207c4762a1bSJed Brown 208c4762a1bSJed Brown PetscFunctionBeginUser; 2099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 210c4762a1bSJed Brown 2119371c9d4SSatish Balay J[0][0] = a; 2129371c9d4SSatish Balay J[0][1] = -1.0; 2139371c9d4SSatish Balay J[1][0] = user->mu * (1.0 + 2.0 * u[0] * u[1]); 2149371c9d4SSatish Balay J[1][1] = a - user->mu * (1.0 - u[0] * u[0]); 2159566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 2169566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 2179566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2189566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 219c4762a1bSJed Brown if (A != B) { 2209566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2219566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 2249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 225*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 226c4762a1bSJed Brown } 227c4762a1bSJed Brown 228d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec U, Mat A, void *ctx) 229d71ae5a4SJacob Faibussowitsch { 230c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {0}; 231c4762a1bSJed Brown PetscScalar J[2][1]; 232c4762a1bSJed Brown const PetscScalar *u; 233c4762a1bSJed Brown 234c4762a1bSJed Brown PetscFunctionBeginUser; 2359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 236c4762a1bSJed Brown 237c4762a1bSJed Brown J[0][0] = 0; 238c4762a1bSJed Brown J[1][0] = (1. - u[0] * u[0]) * u[1] - u[0]; 2399566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES)); 2409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 2419566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2429566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 245*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 246c4762a1bSJed Brown } 247c4762a1bSJed Brown 248d71ae5a4SJacob Faibussowitsch static PetscErrorCode IHessianProductUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 249d71ae5a4SJacob Faibussowitsch { 250c4762a1bSJed Brown const PetscScalar *vl, *vr, *u; 251c4762a1bSJed Brown PetscScalar *vhv; 252c4762a1bSJed Brown PetscScalar dJdU[2][2][2] = {{{0}}}; 253c4762a1bSJed Brown PetscInt i, j, k; 254c4762a1bSJed Brown User user = (User)ctx; 255c4762a1bSJed Brown 256c4762a1bSJed Brown PetscFunctionBeginUser; 2579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 2589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vl[0], &vl)); 2599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vr, &vr)); 2609566063dSJacob Faibussowitsch PetscCall(VecGetArray(VHV[0], &vhv)); 261c4762a1bSJed Brown 262c4762a1bSJed Brown dJdU[1][0][0] = 2. * user->mu * u[1]; 263c4762a1bSJed Brown dJdU[1][1][0] = 2. * user->mu * u[0]; 264c4762a1bSJed Brown dJdU[1][0][1] = 2. * user->mu * u[0]; 265c4762a1bSJed Brown for (j = 0; j < 2; j++) { 266c4762a1bSJed Brown vhv[j] = 0; 267c4762a1bSJed Brown for (k = 0; k < 2; k++) 2689371c9d4SSatish Balay for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k]; 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 2719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 2729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vl[0], &vl)); 2739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vr, &vr)); 2749566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(VHV[0], &vhv)); 275*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 276c4762a1bSJed Brown } 277c4762a1bSJed Brown 278d71ae5a4SJacob Faibussowitsch static PetscErrorCode IHessianProductUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 279d71ae5a4SJacob Faibussowitsch { 280c4762a1bSJed Brown const PetscScalar *vl, *vr, *u; 281c4762a1bSJed Brown PetscScalar *vhv; 282c4762a1bSJed Brown PetscScalar dJdP[2][2][1] = {{{0}}}; 283c4762a1bSJed Brown PetscInt i, j, k; 284c4762a1bSJed Brown 285c4762a1bSJed Brown PetscFunctionBeginUser; 2869566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 2879566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vl[0], &vl)); 2889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vr, &vr)); 2899566063dSJacob Faibussowitsch PetscCall(VecGetArray(VHV[0], &vhv)); 290c4762a1bSJed Brown 291c4762a1bSJed Brown dJdP[1][0][0] = 1. + 2. * u[0] * u[1]; 292c4762a1bSJed Brown dJdP[1][1][0] = u[0] * u[0] - 1.; 293c4762a1bSJed Brown for (j = 0; j < 2; j++) { 294c4762a1bSJed Brown vhv[j] = 0; 295c4762a1bSJed Brown for (k = 0; k < 1; k++) 2969371c9d4SSatish Balay for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdP[i][j][k] * vr[k]; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown 2999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 3009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vl[0], &vl)); 3019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vr, &vr)); 3029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(VHV[0], &vhv)); 303*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 304c4762a1bSJed Brown } 305c4762a1bSJed Brown 306d71ae5a4SJacob Faibussowitsch static PetscErrorCode IHessianProductPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 307d71ae5a4SJacob Faibussowitsch { 308c4762a1bSJed Brown const PetscScalar *vl, *vr, *u; 309c4762a1bSJed Brown PetscScalar *vhv; 310c4762a1bSJed Brown PetscScalar dJdU[2][1][2] = {{{0}}}; 311c4762a1bSJed Brown PetscInt i, j, k; 312c4762a1bSJed Brown 313c4762a1bSJed Brown PetscFunctionBeginUser; 3149566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 3159566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vl[0], &vl)); 3169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Vr, &vr)); 3179566063dSJacob Faibussowitsch PetscCall(VecGetArray(VHV[0], &vhv)); 318c4762a1bSJed Brown 319c4762a1bSJed Brown dJdU[1][0][0] = 1. + 2. * u[1] * u[0]; 320c4762a1bSJed Brown dJdU[1][0][1] = u[0] * u[0] - 1.; 321c4762a1bSJed Brown for (j = 0; j < 1; j++) { 322c4762a1bSJed Brown vhv[j] = 0; 323c4762a1bSJed Brown for (k = 0; k < 2; k++) 3249371c9d4SSatish Balay for (i = 0; i < 2; i++) vhv[j] += vl[i] * dJdU[i][j][k] * vr[k]; 325c4762a1bSJed Brown } 326c4762a1bSJed Brown 3279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 3289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vl[0], &vl)); 3299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Vr, &vr)); 3309566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(VHV[0], &vhv)); 331*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 332c4762a1bSJed Brown } 333c4762a1bSJed Brown 334d71ae5a4SJacob Faibussowitsch static PetscErrorCode IHessianProductPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx) 335d71ae5a4SJacob Faibussowitsch { 336c4762a1bSJed Brown PetscFunctionBeginUser; 337*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 338c4762a1bSJed Brown } 339c4762a1bSJed Brown 340c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 341d71ae5a4SJacob Faibussowitsch static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) 342d71ae5a4SJacob Faibussowitsch { 343c4762a1bSJed Brown const PetscScalar *x; 344c4762a1bSJed Brown PetscReal tfinal, dt; 345c4762a1bSJed Brown User user = (User)ctx; 346c4762a1bSJed Brown Vec interpolatedX; 347c4762a1bSJed Brown 348c4762a1bSJed Brown PetscFunctionBeginUser; 3499566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 3509566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts, &tfinal)); 351c4762a1bSJed Brown 352c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) { 3539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &interpolatedX)); 3549566063dSJacob Faibussowitsch PetscCall(TSInterpolate(ts, user->next_output, interpolatedX)); 3559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(interpolatedX, &x)); 3569371c9d4SSatish Balay PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%g] %" PetscInt_FMT " TS %g (dt = %g) X %g %g\n", (double)user->next_output, step, (double)t, (double)dt, (double)PetscRealPart(x[0]), (double)PetscRealPart(x[1]))); 3579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(interpolatedX, &x)); 3589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&interpolatedX)); 359c4762a1bSJed Brown user->next_output += PetscRealConstant(0.1); 360c4762a1bSJed Brown } 361*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 362c4762a1bSJed Brown } 363c4762a1bSJed Brown 364d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 365d71ae5a4SJacob Faibussowitsch { 366c4762a1bSJed Brown Vec P; 367c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 368c4762a1bSJed Brown PetscScalar *x_ptr; 369c4762a1bSJed Brown const PetscScalar *y_ptr; 370c4762a1bSJed Brown PetscMPIInt size; 371c4762a1bSJed Brown struct _n_User user; 372c4762a1bSJed Brown Tao tao; 373c4762a1bSJed Brown KSP ksp; 374c4762a1bSJed Brown PC pc; 375c4762a1bSJed Brown 376c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 377c4762a1bSJed Brown Initialize program 378c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 379327415f7SBarry Smith PetscFunctionBeginUser; 3809566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 3823c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 383c4762a1bSJed Brown 384c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 3859566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 3869566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBQNLS)); 387c4762a1bSJed Brown 388c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 389c4762a1bSJed Brown Set runtime options 390c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 391c4762a1bSJed Brown user.next_output = 0.0; 392c4762a1bSJed Brown user.mu = PetscRealConstant(1.0e3); 393c4762a1bSJed Brown user.ftime = PetscRealConstant(0.5); 394c4762a1bSJed Brown user.implicitform = PETSC_TRUE; 395c4762a1bSJed Brown 3969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 3979566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL)); 3989566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &user.implicitform, NULL)); 399c4762a1bSJed Brown 400c4762a1bSJed Brown /* Create necessary matrix and vectors, solve same ODE on every process */ 4019566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user.A)); 4029566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 4039566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user.A)); 4049566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.A)); 4059566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &user.U, NULL)); 4069566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &user.Lambda[0], NULL)); 4079566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &user.Lambda2[0], NULL)); 4089566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &user.Ihp1[0], NULL)); 4099566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &user.Ihp2[0], NULL)); 410c4762a1bSJed Brown 4119566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Jacp)); 4129566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1)); 4139566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user.Jacp)); 4149566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.Jacp)); 4159566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.Jacp, &user.Dir, NULL)); 4169566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.Jacp, &user.Mup[0], NULL)); 4179566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.Jacp, &user.Mup2[0], NULL)); 4189566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.Jacp, &user.Ihp3[0], NULL)); 4199566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.Jacp, &user.Ihp4[0], NULL)); 420c4762a1bSJed Brown 421c4762a1bSJed Brown /* Create timestepping solver context */ 4229566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &user.ts)); 4239566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(user.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 424c4762a1bSJed Brown if (user.implicitform) { 4259566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(user.ts, NULL, IFunction, &user)); 4269566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(user.ts, user.A, user.A, IJacobian, &user)); 4279566063dSJacob Faibussowitsch PetscCall(TSSetType(user.ts, TSCN)); 428c4762a1bSJed Brown } else { 4299566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(user.ts, NULL, RHSFunction, &user)); 4309566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(user.ts, user.A, user.A, RHSJacobian, &user)); 4319566063dSJacob Faibussowitsch PetscCall(TSSetType(user.ts, TSRK)); 432c4762a1bSJed Brown } 4339566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(user.ts, user.Jacp, RHSJacobianP, &user)); 4349566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(user.ts, user.ftime)); 4359566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(user.ts, TS_EXACTFINALTIME_MATCHSTEP)); 436c4762a1bSJed Brown 43748a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(user.ts, Monitor, &user, NULL)); 438c4762a1bSJed Brown 439c4762a1bSJed Brown /* Set ODE initial conditions */ 4409566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.U, &x_ptr)); 441c4762a1bSJed Brown x_ptr[0] = 2.0; 442c4762a1bSJed Brown x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu); 4439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user.U, &x_ptr)); 4449566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(user.ts, PetscRealConstant(0.001))); 445c4762a1bSJed Brown 446c4762a1bSJed Brown /* Set runtime options */ 4479566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(user.ts)); 448c4762a1bSJed Brown 4499566063dSJacob Faibussowitsch PetscCall(TSSolve(user.ts, user.U)); 4509566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(user.U, &y_ptr)); 451c4762a1bSJed Brown user.ob[0] = y_ptr[0]; 452c4762a1bSJed Brown user.ob[1] = y_ptr[1]; 4539566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(user.U, &y_ptr)); 454c4762a1bSJed Brown 455c4762a1bSJed Brown /* Save trajectory of solution so that TSAdjointSolve() may be used. 456c4762a1bSJed Brown Skip checkpointing for the first TSSolve since no adjoint run follows it. 457c4762a1bSJed Brown */ 4589566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(user.ts)); 459c4762a1bSJed Brown 460c4762a1bSJed Brown /* Optimization starts */ 4619566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user.H)); 4629566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.H, PETSC_DECIDE, PETSC_DECIDE, 1, 1)); 4639566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.H)); /* Hessian should be symmetric. Do we need to do MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE) ? */ 464c4762a1bSJed Brown 465c4762a1bSJed Brown /* Set initial solution guess */ 4669566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.Jacp, &P, NULL)); 4679566063dSJacob Faibussowitsch PetscCall(VecGetArray(P, &x_ptr)); 468c4762a1bSJed Brown x_ptr[0] = PetscRealConstant(1.2); 4699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(P, &x_ptr)); 4709566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, P)); 471c4762a1bSJed Brown 472c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 4739566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user)); 4749566063dSJacob Faibussowitsch PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user)); 475c4762a1bSJed Brown 476c4762a1bSJed Brown /* Check for any TAO command line options */ 4779566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 478c4762a1bSJed Brown if (ksp) { 4799566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 4809566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 481c4762a1bSJed Brown } 4829566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 483c4762a1bSJed Brown 4849566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 485c4762a1bSJed Brown 4869566063dSJacob Faibussowitsch PetscCall(VecView(P, PETSC_VIEWER_STDOUT_WORLD)); 487c4762a1bSJed Brown /* Free TAO data structures */ 4889566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 489c4762a1bSJed Brown 490c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 491c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 492c4762a1bSJed Brown are no longer needed. 493c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 4949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.H)); 4959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.A)); 4969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.U)); 4979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.Jacp)); 4989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.Lambda[0])); 4999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.Mup[0])); 5009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.Lambda2[0])); 5019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.Mup2[0])); 5029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.Ihp1[0])); 5039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.Ihp2[0])); 5049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.Ihp3[0])); 5059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.Ihp4[0])); 5069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.Dir)); 5079566063dSJacob Faibussowitsch PetscCall(TSDestroy(&user.ts)); 5089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&P)); 5099566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 510b122ec5aSJacob Faibussowitsch return 0; 511c4762a1bSJed Brown } 512c4762a1bSJed Brown 513c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 514c4762a1bSJed Brown /* 515c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient. 516c4762a1bSJed Brown 517c4762a1bSJed Brown Input Parameters: 518c4762a1bSJed Brown tao - the Tao context 519c4762a1bSJed Brown X - the input vector 520a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 521c4762a1bSJed Brown 522c4762a1bSJed Brown Output Parameters: 523c4762a1bSJed Brown f - the newly evaluated function 524c4762a1bSJed Brown G - the newly evaluated gradient 525c4762a1bSJed Brown */ 526d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec P, PetscReal *f, Vec G, void *ctx) 527d71ae5a4SJacob Faibussowitsch { 528c4762a1bSJed Brown User user_ptr = (User)ctx; 529c4762a1bSJed Brown TS ts = user_ptr->ts; 530c4762a1bSJed Brown PetscScalar *x_ptr, *g; 531c4762a1bSJed Brown const PetscScalar *y_ptr; 532c4762a1bSJed Brown 533c4762a1bSJed Brown PetscFunctionBeginUser; 5349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, &y_ptr)); 535c4762a1bSJed Brown user_ptr->mu = y_ptr[0]; 5369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, &y_ptr)); 537c4762a1bSJed Brown 5389566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 5399566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 5409566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, PetscRealConstant(0.001))); /* can be overwritten by command line options */ 5419566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 5429566063dSJacob Faibussowitsch PetscCall(VecGetArray(user_ptr->U, &x_ptr)); 543c4762a1bSJed Brown x_ptr[0] = 2.0; 544c4762a1bSJed Brown x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user_ptr->mu) - 292.0 / (2187.0 * user_ptr->mu * user_ptr->mu); 5459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user_ptr->U, &x_ptr)); 546c4762a1bSJed Brown 5479566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, user_ptr->U)); 548c4762a1bSJed Brown 5499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(user_ptr->U, &y_ptr)); 550c4762a1bSJed Brown *f = (y_ptr[0] - user_ptr->ob[0]) * (y_ptr[0] - user_ptr->ob[0]) + (y_ptr[1] - user_ptr->ob[1]) * (y_ptr[1] - user_ptr->ob[1]); 551c4762a1bSJed Brown 552c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 5539566063dSJacob Faibussowitsch PetscCall(VecGetArray(user_ptr->Lambda[0], &x_ptr)); 554c4762a1bSJed Brown x_ptr[0] = 2. * (y_ptr[0] - user_ptr->ob[0]); 555c4762a1bSJed Brown x_ptr[1] = 2. * (y_ptr[1] - user_ptr->ob[1]); 5569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(user_ptr->U, &y_ptr)); 5579566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user_ptr->Lambda[0], &x_ptr)); 558c4762a1bSJed Brown 5599566063dSJacob Faibussowitsch PetscCall(VecGetArray(user_ptr->Mup[0], &x_ptr)); 560c4762a1bSJed Brown x_ptr[0] = 0.0; 5619566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user_ptr->Mup[0], &x_ptr)); 5629566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, user_ptr->Lambda, user_ptr->Mup)); 563c4762a1bSJed Brown 5649566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 565c4762a1bSJed Brown 5669566063dSJacob Faibussowitsch PetscCall(VecGetArray(user_ptr->Mup[0], &x_ptr)); 5679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(user_ptr->Lambda[0], &y_ptr)); 5689566063dSJacob Faibussowitsch PetscCall(VecGetArray(G, &g)); 569c4762a1bSJed Brown g[0] = y_ptr[1] * (-10.0 / (81.0 * user_ptr->mu * user_ptr->mu) + 2.0 * 292.0 / (2187.0 * user_ptr->mu * user_ptr->mu * user_ptr->mu)) + x_ptr[0]; 5709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user_ptr->Mup[0], &x_ptr)); 5719566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(user_ptr->Lambda[0], &y_ptr)); 5729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G, &g)); 573*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 574c4762a1bSJed Brown } 575c4762a1bSJed Brown 576d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec P, Mat H, Mat Hpre, void *ctx) 577d71ae5a4SJacob Faibussowitsch { 578c4762a1bSJed Brown User user_ptr = (User)ctx; 579c4762a1bSJed Brown PetscScalar harr[1]; 580c4762a1bSJed Brown const PetscInt rows[1] = {0}; 581c4762a1bSJed Brown PetscInt col = 0; 582c4762a1bSJed Brown 583c4762a1bSJed Brown PetscFunctionBeginUser; 5849566063dSJacob Faibussowitsch PetscCall(Adjoint2(P, harr, user_ptr)); 5859566063dSJacob Faibussowitsch PetscCall(MatSetValues(H, 1, rows, 1, &col, harr, INSERT_VALUES)); 586c4762a1bSJed Brown 5879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY)); 5889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY)); 589c4762a1bSJed Brown if (H != Hpre) { 5909566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Hpre, MAT_FINAL_ASSEMBLY)); 5919566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Hpre, MAT_FINAL_ASSEMBLY)); 592c4762a1bSJed Brown } 593*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 594c4762a1bSJed Brown } 595c4762a1bSJed Brown 596d71ae5a4SJacob Faibussowitsch PetscErrorCode Adjoint2(Vec P, PetscScalar arr[], User ctx) 597d71ae5a4SJacob Faibussowitsch { 598c4762a1bSJed Brown TS ts = ctx->ts; 599c4762a1bSJed Brown const PetscScalar *z_ptr; 600c4762a1bSJed Brown PetscScalar *x_ptr, *y_ptr, dzdp, dzdp2; 601c4762a1bSJed Brown Mat tlmsen; 602c4762a1bSJed Brown 603c4762a1bSJed Brown PetscFunctionBeginUser; 604c4762a1bSJed Brown /* Reset TSAdjoint so that AdjointSetUp will be called again */ 6059566063dSJacob Faibussowitsch PetscCall(TSAdjointReset(ts)); 606c4762a1bSJed Brown 607c4762a1bSJed Brown /* The directional vector should be 1 since it is one-dimensional */ 6089566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->Dir, &x_ptr)); 609c4762a1bSJed Brown x_ptr[0] = 1.; 6109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->Dir, &x_ptr)); 611c4762a1bSJed Brown 6129566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, &z_ptr)); 613c4762a1bSJed Brown ctx->mu = z_ptr[0]; 6149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, &z_ptr)); 615c4762a1bSJed Brown 616c4762a1bSJed Brown dzdp = -10.0 / (81.0 * ctx->mu * ctx->mu) + 2.0 * 292.0 / (2187.0 * ctx->mu * ctx->mu * ctx->mu); 617c4762a1bSJed Brown dzdp2 = 2. * 10.0 / (81.0 * ctx->mu * ctx->mu * ctx->mu) - 3.0 * 2.0 * 292.0 / (2187.0 * ctx->mu * ctx->mu * ctx->mu * ctx->mu); 618c4762a1bSJed Brown 6199566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 6209566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 6219566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, PetscRealConstant(0.001))); /* can be overwritten by command line options */ 6229566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 6239566063dSJacob Faibussowitsch PetscCall(TSSetCostHessianProducts(ts, 1, ctx->Lambda2, ctx->Mup2, ctx->Dir)); 624c4762a1bSJed Brown 6259566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(ctx->Jacp)); 6269566063dSJacob Faibussowitsch PetscCall(MatSetValue(ctx->Jacp, 1, 0, dzdp, INSERT_VALUES)); 6279566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ctx->Jacp, MAT_FINAL_ASSEMBLY)); 6289566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ctx->Jacp, MAT_FINAL_ASSEMBLY)); 629c4762a1bSJed Brown 6309566063dSJacob Faibussowitsch PetscCall(TSAdjointSetForward(ts, ctx->Jacp)); 6319566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->U, &y_ptr)); 632c4762a1bSJed Brown y_ptr[0] = 2.0; 633c4762a1bSJed Brown y_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * ctx->mu) - 292.0 / (2187.0 * ctx->mu * ctx->mu); 6349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->U, &y_ptr)); 6359566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, ctx->U)); 636c4762a1bSJed Brown 637c4762a1bSJed Brown /* Set terminal conditions for first- and second-order adjonts */ 6389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ctx->U, &z_ptr)); 6399566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->Lambda[0], &y_ptr)); 640c4762a1bSJed Brown y_ptr[0] = 2. * (z_ptr[0] - ctx->ob[0]); 641c4762a1bSJed Brown y_ptr[1] = 2. * (z_ptr[1] - ctx->ob[1]); 6429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->Lambda[0], &y_ptr)); 6439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ctx->U, &z_ptr)); 6449566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->Mup[0], &y_ptr)); 645c4762a1bSJed Brown y_ptr[0] = 0.0; 6469566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->Mup[0], &y_ptr)); 6479566063dSJacob Faibussowitsch PetscCall(TSForwardGetSensitivities(ts, NULL, &tlmsen)); 6489566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(tlmsen, 0, &x_ptr)); 6499566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->Lambda2[0], &y_ptr)); 650c4762a1bSJed Brown y_ptr[0] = 2. * x_ptr[0]; 651c4762a1bSJed Brown y_ptr[1] = 2. * x_ptr[1]; 6529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->Lambda2[0], &y_ptr)); 6539566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->Mup2[0], &y_ptr)); 654c4762a1bSJed Brown y_ptr[0] = 0.0; 6559566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->Mup2[0], &y_ptr)); 6569566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(tlmsen, &x_ptr)); 6579566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, ctx->Lambda, ctx->Mup)); 658c4762a1bSJed Brown if (ctx->implicitform) { 6599566063dSJacob Faibussowitsch PetscCall(TSSetIHessianProduct(ts, ctx->Ihp1, IHessianProductUU, ctx->Ihp2, IHessianProductUP, ctx->Ihp3, IHessianProductPU, ctx->Ihp4, IHessianProductPP, ctx)); 660c4762a1bSJed Brown } else { 6619566063dSJacob Faibussowitsch PetscCall(TSSetRHSHessianProduct(ts, ctx->Ihp1, RHSHessianProductUU, ctx->Ihp2, RHSHessianProductUP, ctx->Ihp3, RHSHessianProductPU, ctx->Ihp4, RHSHessianProductPP, ctx)); 662c4762a1bSJed Brown } 6639566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 664c4762a1bSJed Brown 6659566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->Lambda[0], &x_ptr)); 6669566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->Lambda2[0], &y_ptr)); 6679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ctx->Mup2[0], &z_ptr)); 668c4762a1bSJed Brown 669c4762a1bSJed Brown arr[0] = x_ptr[1] * dzdp2 + y_ptr[1] * dzdp2 + z_ptr[0]; 670c4762a1bSJed Brown 6719566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->Lambda2[0], &x_ptr)); 6729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->Lambda2[0], &y_ptr)); 6739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ctx->Mup2[0], &z_ptr)); 674c4762a1bSJed Brown 675c4762a1bSJed Brown /* Disable second-order adjoint mode */ 6769566063dSJacob Faibussowitsch PetscCall(TSAdjointReset(ts)); 6779566063dSJacob Faibussowitsch PetscCall(TSAdjointResetForward(ts)); 678*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 679c4762a1bSJed Brown } 680c4762a1bSJed Brown 681c4762a1bSJed Brown /*TEST 682c4762a1bSJed Brown build: 683c4762a1bSJed Brown requires: !complex !single 684c4762a1bSJed Brown test: 685c4762a1bSJed Brown args: -implicitform 0 -ts_type rk -ts_adapt_type none -mu 10 -ts_dt 0.1 -viewer_binary_skip_info -tao_monitor -tao_view 686c4762a1bSJed Brown output_file: output/ex20opt_p_1.out 687c4762a1bSJed Brown 688c4762a1bSJed Brown test: 689c4762a1bSJed Brown suffix: 2 690c4762a1bSJed Brown args: -implicitform 0 -ts_type rk -ts_adapt_type none -mu 10 -ts_dt 0.01 -viewer_binary_skip_info -tao_monitor -tao_type bntr -tao_bnk_pc_type none 691c4762a1bSJed Brown output_file: output/ex20opt_p_2.out 692c4762a1bSJed Brown 693c4762a1bSJed Brown test: 694c4762a1bSJed Brown suffix: 3 695c4762a1bSJed Brown args: -ts_type cn -ts_adapt_type none -mu 100 -ts_dt 0.01 -viewer_binary_skip_info -tao_monitor -tao_view 696c4762a1bSJed Brown output_file: output/ex20opt_p_3.out 697c4762a1bSJed Brown 698c4762a1bSJed Brown test: 699c4762a1bSJed Brown suffix: 4 700c4762a1bSJed Brown args: -ts_type cn -ts_adapt_type none -mu 100 -ts_dt 0.01 -viewer_binary_skip_info -tao_monitor -tao_type bntr -tao_bnk_pc_type none 701c4762a1bSJed Brown output_file: output/ex20opt_p_4.out 702c4762a1bSJed Brown 703c4762a1bSJed Brown TEST*/ 704