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 */ 40c4762a1bSJed Brown static PetscErrorCode RHSFunctionPassive(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 41c4762a1bSJed Brown { 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)); 53c4762a1bSJed Brown PetscFunctionReturn(0); 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 */ 60c4762a1bSJed Brown static PetscErrorCode RHSFunctionActive(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 61c4762a1bSJed Brown { 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 */ 75c4762a1bSJed Brown x_a[0] <<= x[0]; x_a[1] <<= x[1]; /* Mark as independent */ 76c4762a1bSJed Brown f_a[0] = x_a[1]; 77c4762a1bSJed Brown f_a[1] = mu*(1.-x_a[0]*x_a[0])*x_a[1]-x_a[0]; 78c4762a1bSJed Brown f_a[0] >>= f[0]; 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 */ 89c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx) 90c4762a1bSJed Brown { 91c4762a1bSJed Brown User user=(User)ctx; 92a8c08197SHong Zhang const PetscScalar *x; 93c4762a1bSJed Brown 94c4762a1bSJed Brown PetscFunctionBeginUser; 959566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 969566063dSJacob Faibussowitsch PetscCall(PetscAdolcComputeRHSJacobian(1,A,x,user->adctx)); 979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 98c4762a1bSJed Brown PetscFunctionReturn(0); 99c4762a1bSJed Brown } 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* 102c4762a1bSJed Brown Monitor timesteps and use interpolation to output at integer multiples of 0.1 103c4762a1bSJed Brown */ 104c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 105c4762a1bSJed Brown { 106c4762a1bSJed Brown const PetscScalar *x; 107c4762a1bSJed Brown PetscReal tfinal, dt, tprev; 108c4762a1bSJed Brown User user = (User)ctx; 109c4762a1bSJed Brown 110c4762a1bSJed Brown PetscFunctionBeginUser; 1119566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 1129566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts,&tfinal)); 1139566063dSJacob Faibussowitsch PetscCall(TSGetPrevTime(ts,&tprev)); 1149566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 11563a3b9bcSJacob 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]))); 1169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev)); 1179566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 118c4762a1bSJed Brown PetscFunctionReturn(0); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121c4762a1bSJed Brown int main(int argc,char **argv) 122c4762a1bSJed Brown { 123410585f6SHong Zhang TS ts = NULL; /* nonlinear solver */ 124c4762a1bSJed Brown Vec ic,r; 125c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 126c4762a1bSJed Brown PetscScalar *x_ptr; 127c4762a1bSJed Brown PetscMPIInt size; 128c4762a1bSJed Brown struct _n_User user; 129c4762a1bSJed Brown AdolcCtx *adctx; 130c4762a1bSJed Brown Tao tao; 131c4762a1bSJed Brown KSP ksp; 132c4762a1bSJed Brown PC pc; 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135c4762a1bSJed Brown Initialize program 136c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 137*327415f7SBarry Smith PetscFunctionBeginUser; 1389566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 1399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 14008401ef6SPierre Jolivet PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 143c4762a1bSJed Brown Set runtime options and create AdolcCtx 144c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1459566063dSJacob Faibussowitsch PetscCall(PetscNew(&adctx)); 146c4762a1bSJed Brown user.mu = 1.0; 147c4762a1bSJed Brown user.next_output = 0.0; 148c4762a1bSJed Brown user.steps = 0; 149c4762a1bSJed Brown user.ftime = 0.5; 150c4762a1bSJed Brown adctx->m = 2;adctx->n = 2;adctx->p = 2; 151c4762a1bSJed Brown user.adctx = adctx; 152c4762a1bSJed Brown 1539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL)); 1549566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 157c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 158c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1599566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user.A)); 1609566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2)); 1619566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user.A)); 1629566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.A)); 1639566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A,&user.x,NULL)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 166c4762a1bSJed Brown Set initial conditions 167c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1689566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.x,&x_ptr)); 169c4762a1bSJed Brown x_ptr[0] = 2.0; 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)); 190c4762a1bSJed Brown if (monitor) { 1919566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,Monitor,&user,NULL)); 192c4762a1bSJed Brown } 193c4762a1bSJed Brown 1949566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts,0.0)); 19563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %" PetscInt_FMT ", ftime %g\n",(double)user.mu,user.steps,(double)(user.ftime))); 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 198c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used 199c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2009566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 203c4762a1bSJed Brown Set runtime options 204c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2059566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 208c4762a1bSJed Brown Solve nonlinear system 209c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2109566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,user.x)); 2119566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&(user.ftime))); 2129566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&user.steps)); 21363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %" PetscInt_FMT ", ftime %g\n",(double)user.mu,user.steps,(double)user.ftime)); 214c4762a1bSJed Brown 2159566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.x,&x_ptr)); 216c4762a1bSJed Brown user.x_ob[0] = x_ptr[0]; 217c4762a1bSJed Brown user.x_ob[1] = x_ptr[1]; 2189566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user.x,&x_ptr)); 219c4762a1bSJed Brown 2209566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A,&user.lambda[0],NULL)); 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 2239566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao)); 2249566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao,TAOCG)); 225c4762a1bSJed Brown 226c4762a1bSJed Brown /* Set initial solution guess */ 2279566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A,&ic,NULL)); 2289566063dSJacob Faibussowitsch PetscCall(VecGetArray(ic,&x_ptr)); 229c4762a1bSJed Brown x_ptr[0] = 2.1; 230c4762a1bSJed Brown x_ptr[1] = 0.7; 2319566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ic,&x_ptr)); 232c4762a1bSJed Brown 2339566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao,ic)); 234c4762a1bSJed Brown 235c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 2369566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user)); 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* Check for any TAO command line options */ 2399566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 2409566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao,&ksp)); 241c4762a1bSJed Brown if (ksp) { 2429566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 2439566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCNONE)); 244c4762a1bSJed Brown } 245c4762a1bSJed Brown 2469566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(tao,1e-10,PETSC_DEFAULT,PETSC_DEFAULT)); 247c4762a1bSJed Brown 248c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 2499566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 250c4762a1bSJed Brown 251c4762a1bSJed Brown /* Free TAO data structures */ 2529566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 255c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 256c4762a1bSJed Brown are no longer needed. 257c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.A)); 2599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.x)); 2609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.lambda[0])); 2619566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ic)); 2639566063dSJacob Faibussowitsch PetscCall(PetscFree(adctx)); 2649566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 265b122ec5aSJacob Faibussowitsch return 0; 266c4762a1bSJed Brown } 267c4762a1bSJed Brown 268c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 269c4762a1bSJed Brown /* 270c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient. 271c4762a1bSJed Brown 272c4762a1bSJed Brown Input Parameters: 273c4762a1bSJed Brown tao - the Tao context 274c4762a1bSJed Brown X - the input vector 275a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 276c4762a1bSJed Brown 277c4762a1bSJed Brown Output Parameters: 278c4762a1bSJed Brown f - the newly evaluated function 279c4762a1bSJed Brown G - the newly evaluated gradient 280c4762a1bSJed Brown */ 281c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx) 282c4762a1bSJed Brown { 283c4762a1bSJed Brown User user = (User)ctx; 284c4762a1bSJed Brown TS ts; 285c4762a1bSJed Brown PetscScalar *x_ptr,*y_ptr; 286c4762a1bSJed Brown 287c4762a1bSJed Brown PetscFunctionBeginUser; 2889566063dSJacob Faibussowitsch PetscCall(VecCopy(IC,user->x)); 289c4762a1bSJed Brown 290c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 291c4762a1bSJed Brown Create timestepping solver context 292c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2939566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 2949566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSRK)); 2959566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,user)); 296c4762a1bSJed Brown /* Set RHS Jacobian for the adjoint integration */ 2979566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts,user->A,user->A,RHSJacobian,user)); 298c4762a1bSJed Brown 299c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 300c4762a1bSJed Brown Set time 301c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3029566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts,0.0)); 3039566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,.001)); 3049566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,0.5)); 3059566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 306c4762a1bSJed Brown 3079566063dSJacob Faibussowitsch PetscCall(TSSetTolerances(ts,1e-7,NULL,1e-7,NULL)); 308c4762a1bSJed Brown 309c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 310c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used 311c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3129566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 313c4762a1bSJed Brown 314c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 315c4762a1bSJed Brown Set runtime options 316c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3179566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 318c4762a1bSJed Brown 3199566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,user->x)); 3209566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&user->ftime)); 3219566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&user->steps)); 32263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %.6f, steps %" PetscInt_FMT ", ftime %g\n",(double)user->mu,user->steps,(double)user->ftime)); 323c4762a1bSJed Brown 3249566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->x,&x_ptr)); 325c4762a1bSJed 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]); 3269566063dSJacob 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))); 327c4762a1bSJed Brown 328c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 329c4762a1bSJed Brown Adjoint model starts here 330c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 331c4762a1bSJed Brown /* Redet initial conditions for the adjoint integration */ 3329566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->lambda[0],&y_ptr)); 333c4762a1bSJed Brown y_ptr[0] = 2.*(x_ptr[0]-user->x_ob[0]); 334c4762a1bSJed Brown y_ptr[1] = 2.*(x_ptr[1]-user->x_ob[1]); 3359566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->lambda[0],&y_ptr)); 3369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->x,&x_ptr)); 3379566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts,1,user->lambda,NULL)); 338c4762a1bSJed Brown 3399566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 340c4762a1bSJed Brown 3419566063dSJacob Faibussowitsch PetscCall(VecCopy(user->lambda[0],G)); 342c4762a1bSJed Brown 3439566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 344c4762a1bSJed Brown PetscFunctionReturn(0); 345c4762a1bSJed Brown } 346c4762a1bSJed Brown 347c4762a1bSJed Brown /*TEST 348c4762a1bSJed Brown 349c4762a1bSJed Brown build: 350c4762a1bSJed Brown requires: double !complex adolc 351c4762a1bSJed Brown 352c4762a1bSJed Brown test: 353c4762a1bSJed Brown suffix: 1 354c4762a1bSJed Brown args: -ts_rhs_jacobian_test_mult_transpose FALSE -tao_max_it 2 -ts_rhs_jacobian_test_mult FALSE 355c4762a1bSJed Brown output_file: output/ex16opt_ic_1.out 356c4762a1bSJed Brown 357c4762a1bSJed Brown TEST*/ 358