1c4762a1bSJed Brown static char help[] = "Demonstrates automatic Jacobian generation using ADOL-C for an ODE-constrained optimization problem.\n\ 2c4762a1bSJed Brown Input parameters include:\n\ 3c4762a1bSJed Brown -mu : stiffness parameter\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown /* 6c4762a1bSJed Brown REQUIRES configuration of PETSc with option --download-adolc. 7c4762a1bSJed Brown 8c4762a1bSJed Brown For documentation on ADOL-C, see 9c4762a1bSJed Brown $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown /* ------------------------------------------------------------------------ 12c4762a1bSJed Brown See ex16opt_ic for a description of the problem being solved. 13c4762a1bSJed Brown ------------------------------------------------------------------------- */ 14c4762a1bSJed Brown #include <petsctao.h> 15c4762a1bSJed Brown #include <petscts.h> 16c4762a1bSJed Brown #include <petscmat.h> 17c4762a1bSJed Brown #include "adolc-utils/drivers.cxx" 18c4762a1bSJed Brown #include <adolc/adolc.h> 19c4762a1bSJed Brown 20c4762a1bSJed Brown typedef struct _n_User *User; 21c4762a1bSJed Brown struct _n_User { 22c4762a1bSJed Brown PetscReal mu; 23c4762a1bSJed Brown PetscReal next_output; 24c4762a1bSJed Brown PetscInt steps; 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* Sensitivity analysis support */ 27c4762a1bSJed Brown PetscReal ftime, x_ob[2]; 28c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 29c4762a1bSJed Brown Vec x, lambda[2]; /* adjoint variables */ 30c4762a1bSJed Brown 31c4762a1bSJed Brown /* Automatic differentiation support */ 32c4762a1bSJed Brown AdolcCtx *adctx; 33c4762a1bSJed Brown }; 34c4762a1bSJed Brown 35c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* 38c4762a1bSJed Brown 'Passive' RHS function, used in residual evaluations during the time integration. 39c4762a1bSJed Brown */ 409371c9d4SSatish Balay static PetscErrorCode RHSFunctionPassive(TS ts, PetscReal t, Vec X, Vec F, void *ctx) { 41c4762a1bSJed Brown User user = (User)ctx; 42c4762a1bSJed Brown PetscScalar *f; 43c4762a1bSJed Brown const PetscScalar *x; 44c4762a1bSJed Brown 45c4762a1bSJed Brown PetscFunctionBeginUser; 469566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 479566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 48c4762a1bSJed Brown f[0] = x[1]; 49c4762a1bSJed Brown f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0]; 509566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 52c4762a1bSJed Brown PetscFunctionReturn(0); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* 56c4762a1bSJed Brown Trace RHS to mark on tape 1 the dependence of f upon x. This tape is used in generating the 57c4762a1bSJed Brown Jacobian transform. 58c4762a1bSJed Brown */ 599371c9d4SSatish Balay static PetscErrorCode RHSFunctionActive(TS ts, PetscReal t, Vec X, Vec F, void *ctx) { 60c4762a1bSJed Brown User user = (User)ctx; 61c4762a1bSJed Brown PetscReal mu = user->mu; 62c4762a1bSJed Brown PetscScalar *f; 63c4762a1bSJed Brown const PetscScalar *x; 64c4762a1bSJed Brown 65c4762a1bSJed Brown adouble f_a[2]; /* adouble for dependent variables */ 66c4762a1bSJed Brown adouble x_a[2]; /* adouble for independent variables */ 67c4762a1bSJed Brown 68c4762a1bSJed Brown PetscFunctionBeginUser; 699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 709566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 71c4762a1bSJed Brown 72c4762a1bSJed Brown trace_on(1); /* Start of active section */ 739371c9d4SSatish Balay x_a[0] <<= x[0]; 749371c9d4SSatish Balay x_a[1] <<= x[1]; /* Mark as independent */ 75c4762a1bSJed Brown f_a[0] = x_a[1]; 76c4762a1bSJed Brown f_a[1] = mu * (1. - x_a[0] * x_a[0]) * x_a[1] - x_a[0]; 779371c9d4SSatish Balay f_a[0] >>= f[0]; 789371c9d4SSatish Balay f_a[1] >>= f[1]; /* Mark as dependent */ 79c4762a1bSJed Brown trace_off(1); /* End of active section */ 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 83c4762a1bSJed Brown PetscFunctionReturn(0); 84c4762a1bSJed Brown } 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* 87c4762a1bSJed Brown Compute the Jacobian w.r.t. x using PETSc-ADOL-C driver. 88c4762a1bSJed Brown */ 899371c9d4SSatish Balay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat A, Mat B, void *ctx) { 90c4762a1bSJed Brown User user = (User)ctx; 91a8c08197SHong Zhang const PetscScalar *x; 92c4762a1bSJed Brown 93c4762a1bSJed Brown PetscFunctionBeginUser; 949566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 959566063dSJacob Faibussowitsch PetscCall(PetscAdolcComputeRHSJacobian(1, A, x, user->adctx)); 969566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 97c4762a1bSJed Brown PetscFunctionReturn(0); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* 101c4762a1bSJed Brown Monitor timesteps and use interpolation to output at integer multiples of 0.1 102c4762a1bSJed Brown */ 1039371c9d4SSatish Balay static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) { 104c4762a1bSJed Brown const PetscScalar *x; 105c4762a1bSJed Brown PetscReal tfinal, dt, tprev; 106c4762a1bSJed Brown User user = (User)ctx; 107c4762a1bSJed Brown 108c4762a1bSJed Brown PetscFunctionBeginUser; 1099566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 1109566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts, &tfinal)); 1119566063dSJacob Faibussowitsch PetscCall(TSGetPrevTime(ts, &tprev)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 11363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%.1f] %" PetscInt_FMT " TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n", (double)user->next_output, step, (double)t, (double)dt, (double)PetscRealPart(x[0]), (double)PetscRealPart(x[1]))); 1149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "t %.6f (tprev = %.6f) \n", (double)t, (double)tprev)); 1159566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 116c4762a1bSJed Brown PetscFunctionReturn(0); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 1199371c9d4SSatish Balay int main(int argc, char **argv) { 120410585f6SHong Zhang TS ts = NULL; /* nonlinear solver */ 121c4762a1bSJed Brown Vec ic, r; 122c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 123c4762a1bSJed Brown PetscScalar *x_ptr; 124c4762a1bSJed Brown PetscMPIInt size; 125c4762a1bSJed Brown struct _n_User user; 126c4762a1bSJed Brown AdolcCtx *adctx; 127c4762a1bSJed Brown Tao tao; 128c4762a1bSJed Brown KSP ksp; 129c4762a1bSJed Brown PC pc; 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 132c4762a1bSJed Brown Initialize program 133c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 134327415f7SBarry Smith PetscFunctionBeginUser; 1359566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 13708401ef6SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 140c4762a1bSJed Brown Set runtime options and create AdolcCtx 141c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1429566063dSJacob Faibussowitsch PetscCall(PetscNew(&adctx)); 143c4762a1bSJed Brown user.mu = 1.0; 144c4762a1bSJed Brown user.next_output = 0.0; 145c4762a1bSJed Brown user.steps = 0; 146c4762a1bSJed Brown user.ftime = 0.5; 1479371c9d4SSatish Balay adctx->m = 2; 1489371c9d4SSatish Balay adctx->n = 2; 1499371c9d4SSatish Balay adctx->p = 2; 150c4762a1bSJed Brown user.adctx = adctx; 151c4762a1bSJed Brown 1529566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL)); 1539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 157c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1589566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user.A)); 1599566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 1609566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user.A)); 1619566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.A)); 1629566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &user.x, NULL)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 165c4762a1bSJed Brown Set initial conditions 166c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1679566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.x, &x_ptr)); 1689371c9d4SSatish Balay x_ptr[0] = 2.0; 1699371c9d4SSatish Balay x_ptr[1] = 0.66666654321; 1709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user.x, &x_ptr)); 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 173c4762a1bSJed Brown Trace just once on each tape and put zeros on Jacobian diagonal 174c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1759566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.x, &r)); 1769566063dSJacob Faibussowitsch PetscCall(RHSFunctionActive(ts, 0., user.x, r, &user)); 1779566063dSJacob Faibussowitsch PetscCall(VecSet(r, 0)); 1789566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(user.A, r, INSERT_VALUES)); 1799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 182c4762a1bSJed Brown Create timestepping solver context 183c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1849566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1859566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK)); 1869566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionPassive, &user)); 1879566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, user.A, user.A, RHSJacobian, &user)); 1889566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, user.ftime)); 1899566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 190*48a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); 191c4762a1bSJed Brown 1929566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 19363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, user.steps, (double)(user.ftime))); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 196c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used 197c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1989566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 199c4762a1bSJed Brown 200c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 201c4762a1bSJed Brown Set runtime options 202c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2039566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 206c4762a1bSJed Brown Solve nonlinear system 207c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2089566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, user.x)); 2099566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &(user.ftime))); 2109566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &user.steps)); 21163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, user.steps, (double)user.ftime)); 212c4762a1bSJed Brown 2139566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.x, &x_ptr)); 214c4762a1bSJed Brown user.x_ob[0] = x_ptr[0]; 215c4762a1bSJed Brown user.x_ob[1] = x_ptr[1]; 2169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user.x, &x_ptr)); 217c4762a1bSJed Brown 2189566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &user.lambda[0], NULL)); 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 2219566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 2229566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOCG)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* Set initial solution guess */ 2259566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &ic, NULL)); 2269566063dSJacob Faibussowitsch PetscCall(VecGetArray(ic, &x_ptr)); 227c4762a1bSJed Brown x_ptr[0] = 2.1; 228c4762a1bSJed Brown x_ptr[1] = 0.7; 2299566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ic, &x_ptr)); 230c4762a1bSJed Brown 2319566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, ic)); 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 2349566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user)); 235c4762a1bSJed Brown 236c4762a1bSJed Brown /* Check for any TAO command line options */ 2379566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 2389566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 239c4762a1bSJed Brown if (ksp) { 2409566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 2419566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 242c4762a1bSJed Brown } 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(tao, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT)); 245c4762a1bSJed Brown 246c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 2479566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 248c4762a1bSJed Brown 249c4762a1bSJed Brown /* Free TAO data structures */ 2509566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 251c4762a1bSJed Brown 252c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 253c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 254c4762a1bSJed Brown are no longer needed. 255c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.A)); 2579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.x)); 2589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.lambda[0])); 2599566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ic)); 2619566063dSJacob Faibussowitsch PetscCall(PetscFree(adctx)); 2629566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 263b122ec5aSJacob Faibussowitsch return 0; 264c4762a1bSJed Brown } 265c4762a1bSJed Brown 266c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 267c4762a1bSJed Brown /* 268c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient. 269c4762a1bSJed Brown 270c4762a1bSJed Brown Input Parameters: 271c4762a1bSJed Brown tao - the Tao context 272c4762a1bSJed Brown X - the input vector 273a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 274c4762a1bSJed Brown 275c4762a1bSJed Brown Output Parameters: 276c4762a1bSJed Brown f - the newly evaluated function 277c4762a1bSJed Brown G - the newly evaluated gradient 278c4762a1bSJed Brown */ 2799371c9d4SSatish Balay PetscErrorCode FormFunctionGradient(Tao tao, Vec IC, PetscReal *f, Vec G, void *ctx) { 280c4762a1bSJed Brown User user = (User)ctx; 281c4762a1bSJed Brown TS ts; 282c4762a1bSJed Brown PetscScalar *x_ptr, *y_ptr; 283c4762a1bSJed Brown 284c4762a1bSJed Brown PetscFunctionBeginUser; 2859566063dSJacob Faibussowitsch PetscCall(VecCopy(IC, user->x)); 286c4762a1bSJed Brown 287c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 288c4762a1bSJed Brown Create timestepping solver context 289c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2909566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 2919566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK)); 2929566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionPassive, user)); 293c4762a1bSJed Brown /* Set RHS Jacobian for the adjoint integration */ 2949566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, user->A, user->A, RHSJacobian, user)); 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 297c4762a1bSJed Brown Set time 298c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2999566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 3009566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .001)); 3019566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 0.5)); 3029566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 303c4762a1bSJed Brown 3049566063dSJacob Faibussowitsch PetscCall(TSSetTolerances(ts, 1e-7, NULL, 1e-7, NULL)); 305c4762a1bSJed Brown 306c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 307c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used 308c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3099566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 310c4762a1bSJed Brown 311c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 312c4762a1bSJed Brown Set runtime options 313c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3149566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 315c4762a1bSJed Brown 3169566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, user->x)); 3179566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &user->ftime)); 3189566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &user->steps)); 31963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %.6f, steps %" PetscInt_FMT ", ftime %g\n", (double)user->mu, user->steps, (double)user->ftime)); 320c4762a1bSJed Brown 3219566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->x, &x_ptr)); 322c4762a1bSJed Brown *f = (x_ptr[0] - user->x_ob[0]) * (x_ptr[0] - user->x_ob[0]) + (x_ptr[1] - user->x_ob[1]) * (x_ptr[1] - user->x_ob[1]); 3239566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Observed value y_ob=[%f; %f], ODE solution y=[%f;%f], Cost function f=%f\n", (double)user->x_ob[0], (double)user->x_ob[1], (double)x_ptr[0], (double)x_ptr[1], (double)(*f))); 324c4762a1bSJed Brown 325c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 326c4762a1bSJed Brown Adjoint model starts here 327c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 328c4762a1bSJed Brown /* Redet initial conditions for the adjoint integration */ 3299566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->lambda[0], &y_ptr)); 330c4762a1bSJed Brown y_ptr[0] = 2. * (x_ptr[0] - user->x_ob[0]); 331c4762a1bSJed Brown y_ptr[1] = 2. * (x_ptr[1] - user->x_ob[1]); 3329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->lambda[0], &y_ptr)); 3339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->x, &x_ptr)); 3349566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, user->lambda, NULL)); 335c4762a1bSJed Brown 3369566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 337c4762a1bSJed Brown 3389566063dSJacob Faibussowitsch PetscCall(VecCopy(user->lambda[0], G)); 339c4762a1bSJed Brown 3409566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 341c4762a1bSJed Brown PetscFunctionReturn(0); 342c4762a1bSJed Brown } 343c4762a1bSJed Brown 344c4762a1bSJed Brown /*TEST 345c4762a1bSJed Brown 346c4762a1bSJed Brown build: 347c4762a1bSJed Brown requires: double !complex adolc 348c4762a1bSJed Brown 349c4762a1bSJed Brown test: 350c4762a1bSJed Brown suffix: 1 351c4762a1bSJed Brown args: -ts_rhs_jacobian_test_mult_transpose FALSE -tao_max_it 2 -ts_rhs_jacobian_test_mult FALSE 352c4762a1bSJed Brown output_file: output/ex16opt_ic_1.out 353c4762a1bSJed Brown 354c4762a1bSJed Brown TEST*/ 355