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 */ 40d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunctionPassive(TS ts, PetscReal t, Vec X, Vec F, void *ctx) 41d71ae5a4SJacob Faibussowitsch { 42c4762a1bSJed Brown User user = (User)ctx; 43c4762a1bSJed Brown PetscScalar *f; 44c4762a1bSJed Brown const PetscScalar *x; 45c4762a1bSJed Brown 46c4762a1bSJed Brown PetscFunctionBeginUser; 479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 489566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 49c4762a1bSJed Brown f[0] = x[1]; 50c4762a1bSJed Brown f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0]; 519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* 57c4762a1bSJed Brown Trace RHS to mark on tape 1 the dependence of f upon x. This tape is used in generating the 58c4762a1bSJed Brown Jacobian transform. 59c4762a1bSJed Brown */ 60d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunctionActive(TS ts, PetscReal t, Vec X, Vec F, void *ctx) 61d71ae5a4SJacob Faibussowitsch { 62c4762a1bSJed Brown User user = (User)ctx; 63c4762a1bSJed Brown PetscReal mu = user->mu; 64c4762a1bSJed Brown PetscScalar *f; 65c4762a1bSJed Brown const PetscScalar *x; 66c4762a1bSJed Brown 67c4762a1bSJed Brown adouble f_a[2]; /* adouble for dependent variables */ 68c4762a1bSJed Brown adouble x_a[2]; /* adouble for independent variables */ 69c4762a1bSJed Brown 70c4762a1bSJed Brown PetscFunctionBeginUser; 719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 729566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 73c4762a1bSJed Brown 74c4762a1bSJed Brown trace_on(1); /* Start of active section */ 759371c9d4SSatish Balay x_a[0] <<= x[0]; 769371c9d4SSatish Balay x_a[1] <<= x[1]; /* Mark as independent */ 77c4762a1bSJed Brown f_a[0] = x_a[1]; 78c4762a1bSJed Brown f_a[1] = mu * (1. - x_a[0] * x_a[0]) * x_a[1] - x_a[0]; 799371c9d4SSatish Balay f_a[0] >>= f[0]; 809371c9d4SSatish Balay f_a[1] >>= f[1]; /* Mark as dependent */ 81c4762a1bSJed Brown trace_off(1); /* End of active section */ 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88c4762a1bSJed Brown /* 89c4762a1bSJed Brown Compute the Jacobian w.r.t. x using PETSc-ADOL-C driver. 90c4762a1bSJed Brown */ 91d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat A, Mat B, void *ctx) 92d71ae5a4SJacob Faibussowitsch { 93c4762a1bSJed Brown User user = (User)ctx; 94a8c08197SHong Zhang const PetscScalar *x; 95c4762a1bSJed Brown 96c4762a1bSJed Brown PetscFunctionBeginUser; 979566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 989566063dSJacob Faibussowitsch PetscCall(PetscAdolcComputeRHSJacobian(1, A, x, user->adctx)); 999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* 104c4762a1bSJed Brown Monitor timesteps and use interpolation to output at integer multiples of 0.1 105c4762a1bSJed Brown */ 106d71ae5a4SJacob Faibussowitsch static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) 107d71ae5a4SJacob Faibussowitsch { 108c4762a1bSJed Brown const PetscScalar *x; 109c4762a1bSJed Brown PetscReal tfinal, dt, tprev; 110c4762a1bSJed Brown User user = (User)ctx; 111c4762a1bSJed Brown 112c4762a1bSJed Brown PetscFunctionBeginUser; 1139566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 1149566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts, &tfinal)); 1159566063dSJacob Faibussowitsch PetscCall(TSGetPrevTime(ts, &tprev)); 1169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 11763a3b9bcSJacob 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]))); 1189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "t %.6f (tprev = %.6f) \n", (double)t, (double)tprev)); 1199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 1203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 124d71ae5a4SJacob Faibussowitsch { 125410585f6SHong Zhang TS ts = NULL; /* nonlinear solver */ 126c4762a1bSJed Brown Vec ic, r; 127c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 128c4762a1bSJed Brown PetscScalar *x_ptr; 129c4762a1bSJed Brown PetscMPIInt size; 130c4762a1bSJed Brown struct _n_User user; 131c4762a1bSJed Brown AdolcCtx *adctx; 132c4762a1bSJed Brown Tao tao; 133c4762a1bSJed Brown KSP ksp; 134c4762a1bSJed Brown PC pc; 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 137c4762a1bSJed Brown Initialize program 138c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 139327415f7SBarry Smith PetscFunctionBeginUser; 1409566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 14208401ef6SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 145c4762a1bSJed Brown Set runtime options and create AdolcCtx 146c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1479566063dSJacob Faibussowitsch PetscCall(PetscNew(&adctx)); 148c4762a1bSJed Brown user.mu = 1.0; 149c4762a1bSJed Brown user.next_output = 0.0; 150c4762a1bSJed Brown user.steps = 0; 151c4762a1bSJed Brown user.ftime = 0.5; 1529371c9d4SSatish Balay adctx->m = 2; 1539371c9d4SSatish Balay adctx->n = 2; 1549371c9d4SSatish Balay adctx->p = 2; 155c4762a1bSJed Brown user.adctx = adctx; 156c4762a1bSJed Brown 1579566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL)); 1589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 159c4762a1bSJed Brown 160c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 161c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 162c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1639566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user.A)); 1649566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 1659566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user.A)); 1669566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.A)); 1679566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &user.x, NULL)); 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 170c4762a1bSJed Brown Set initial conditions 171c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1729566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.x, &x_ptr)); 1739371c9d4SSatish Balay x_ptr[0] = 2.0; 1749371c9d4SSatish Balay x_ptr[1] = 0.66666654321; 1759566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user.x, &x_ptr)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 178c4762a1bSJed Brown Trace just once on each tape and put zeros on Jacobian diagonal 179c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(user.x, &r)); 1819566063dSJacob Faibussowitsch PetscCall(RHSFunctionActive(ts, 0., user.x, r, &user)); 1829566063dSJacob Faibussowitsch PetscCall(VecSet(r, 0)); 1839566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(user.A, r, INSERT_VALUES)); 1849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 187c4762a1bSJed Brown Create timestepping solver context 188c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1899566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1909566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK)); 1919566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionPassive, &user)); 1929566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, user.A, user.A, RHSJacobian, &user)); 1939566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, user.ftime)); 1949566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 19548a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); 196c4762a1bSJed Brown 1979566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 198*f4f49eeaSPierre Jolivet PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, user.steps, (double)user.ftime)); 199c4762a1bSJed Brown 200c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 201c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used 202c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2039566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 206c4762a1bSJed Brown Set runtime options 207c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2089566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 211c4762a1bSJed Brown Solve nonlinear system 212c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2139566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, user.x)); 214*f4f49eeaSPierre Jolivet PetscCall(TSGetSolveTime(ts, &user.ftime)); 2159566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &user.steps)); 21663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, user.steps, (double)user.ftime)); 217c4762a1bSJed Brown 2189566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.x, &x_ptr)); 219c4762a1bSJed Brown user.x_ob[0] = x_ptr[0]; 220c4762a1bSJed Brown user.x_ob[1] = x_ptr[1]; 2219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user.x, &x_ptr)); 222c4762a1bSJed Brown 2239566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &user.lambda[0], NULL)); 224c4762a1bSJed Brown 225c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 2269566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 2279566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOCG)); 228c4762a1bSJed Brown 229c4762a1bSJed Brown /* Set initial solution guess */ 2309566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &ic, NULL)); 2319566063dSJacob Faibussowitsch PetscCall(VecGetArray(ic, &x_ptr)); 232c4762a1bSJed Brown x_ptr[0] = 2.1; 233c4762a1bSJed Brown x_ptr[1] = 0.7; 2349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ic, &x_ptr)); 235c4762a1bSJed Brown 2369566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, ic)); 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 2399566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user)); 240c4762a1bSJed Brown 241c4762a1bSJed Brown /* Check for any TAO command line options */ 2429566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 2439566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 244c4762a1bSJed Brown if (ksp) { 2459566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 2469566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 247c4762a1bSJed Brown } 248c4762a1bSJed Brown 2499566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(tao, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT)); 250c4762a1bSJed Brown 251c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 2529566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* Free TAO data structures */ 2559566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 256c4762a1bSJed Brown 257c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 258c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 259c4762a1bSJed Brown are no longer needed. 260c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.A)); 2629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.x)); 2639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.lambda[0])); 2649566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ic)); 2669566063dSJacob Faibussowitsch PetscCall(PetscFree(adctx)); 2679566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 268b122ec5aSJacob Faibussowitsch return 0; 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 271c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 272c4762a1bSJed Brown /* 273c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient. 274c4762a1bSJed Brown 275c4762a1bSJed Brown Input Parameters: 276c4762a1bSJed Brown tao - the Tao context 277c4762a1bSJed Brown X - the input vector 278a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 279c4762a1bSJed Brown 280c4762a1bSJed Brown Output Parameters: 281c4762a1bSJed Brown f - the newly evaluated function 282c4762a1bSJed Brown G - the newly evaluated gradient 283c4762a1bSJed Brown */ 284d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec IC, PetscReal *f, Vec G, void *ctx) 285d71ae5a4SJacob Faibussowitsch { 286c4762a1bSJed Brown User user = (User)ctx; 287c4762a1bSJed Brown TS ts; 288c4762a1bSJed Brown PetscScalar *x_ptr, *y_ptr; 289c4762a1bSJed Brown 290c4762a1bSJed Brown PetscFunctionBeginUser; 2919566063dSJacob Faibussowitsch PetscCall(VecCopy(IC, user->x)); 292c4762a1bSJed Brown 293c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 294c4762a1bSJed Brown Create timestepping solver context 295c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2969566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 2979566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK)); 2989566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionPassive, user)); 299c4762a1bSJed Brown /* Set RHS Jacobian for the adjoint integration */ 3009566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, user->A, user->A, RHSJacobian, user)); 301c4762a1bSJed Brown 302c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 303c4762a1bSJed Brown Set time 304c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3059566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 3069566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .001)); 3079566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 0.5)); 3089566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 309c4762a1bSJed Brown 3109566063dSJacob Faibussowitsch PetscCall(TSSetTolerances(ts, 1e-7, NULL, 1e-7, NULL)); 311c4762a1bSJed Brown 312c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 313c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used 314c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3159566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 316c4762a1bSJed Brown 317c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 318c4762a1bSJed Brown Set runtime options 319c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3209566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 321c4762a1bSJed Brown 3229566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, user->x)); 3239566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &user->ftime)); 3249566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &user->steps)); 32563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %.6f, steps %" PetscInt_FMT ", ftime %g\n", (double)user->mu, user->steps, (double)user->ftime)); 326c4762a1bSJed Brown 3279566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->x, &x_ptr)); 328c4762a1bSJed 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]); 3299566063dSJacob 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))); 330c4762a1bSJed Brown 331c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 332c4762a1bSJed Brown Adjoint model starts here 333c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 334c4762a1bSJed Brown /* Redet initial conditions for the adjoint integration */ 3359566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->lambda[0], &y_ptr)); 336c4762a1bSJed Brown y_ptr[0] = 2. * (x_ptr[0] - user->x_ob[0]); 337c4762a1bSJed Brown y_ptr[1] = 2. * (x_ptr[1] - user->x_ob[1]); 3389566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->lambda[0], &y_ptr)); 3399566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->x, &x_ptr)); 3409566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, user->lambda, NULL)); 341c4762a1bSJed Brown 3429566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 343c4762a1bSJed Brown 3449566063dSJacob Faibussowitsch PetscCall(VecCopy(user->lambda[0], G)); 345c4762a1bSJed Brown 3469566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 3473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown 350c4762a1bSJed Brown /*TEST 351c4762a1bSJed Brown 352c4762a1bSJed Brown build: 353c4762a1bSJed Brown requires: double !complex adolc 354c4762a1bSJed Brown 355c4762a1bSJed Brown test: 356c4762a1bSJed Brown suffix: 1 357c4762a1bSJed Brown args: -ts_rhs_jacobian_test_mult_transpose FALSE -tao_max_it 2 -ts_rhs_jacobian_test_mult FALSE 358c4762a1bSJed Brown output_file: output/ex16opt_ic_1.out 359c4762a1bSJed Brown 360c4762a1bSJed Brown TEST*/ 361