1 static char help[] = "Demonstrates automatic Jacobian generation using ADOL-C for an ODE-constrained optimization problem.\n\ 2 Input parameters include:\n\ 3 -mu : stiffness parameter\n\n"; 4 5 /* 6 REQUIRES configuration of PETSc with option --download-adolc. 7 8 For documentation on ADOL-C, see 9 $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf 10 */ 11 /* ------------------------------------------------------------------------ 12 See ex16opt_ic for a description of the problem being solved. 13 ------------------------------------------------------------------------- */ 14 #include <petsctao.h> 15 #include <petscts.h> 16 #include <petscmat.h> 17 #include "adolc-utils/drivers.cxx" 18 #include <adolc/adolc.h> 19 20 typedef struct _n_User *User; 21 struct _n_User { 22 PetscReal mu; 23 PetscReal next_output; 24 PetscInt steps; 25 26 /* Sensitivity analysis support */ 27 PetscReal ftime,x_ob[2]; 28 Mat A; /* Jacobian matrix */ 29 Vec x,lambda[2]; /* adjoint variables */ 30 31 /* Automatic differentiation support */ 32 AdolcCtx *adctx; 33 }; 34 35 PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*); 36 37 /* 38 'Passive' RHS function, used in residual evaluations during the time integration. 39 */ 40 static PetscErrorCode RHSFunctionPassive(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 41 { 42 User user = (User)ctx; 43 PetscScalar *f; 44 const PetscScalar *x; 45 46 PetscFunctionBeginUser; 47 PetscCall(VecGetArrayRead(X,&x)); 48 PetscCall(VecGetArray(F,&f)); 49 f[0] = x[1]; 50 f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0]; 51 PetscCall(VecRestoreArrayRead(X,&x)); 52 PetscCall(VecRestoreArray(F,&f)); 53 PetscFunctionReturn(0); 54 } 55 56 /* 57 Trace RHS to mark on tape 1 the dependence of f upon x. This tape is used in generating the 58 Jacobian transform. 59 */ 60 static PetscErrorCode RHSFunctionActive(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 61 { 62 User user = (User)ctx; 63 PetscReal mu = user->mu; 64 PetscScalar *f; 65 const PetscScalar *x; 66 67 adouble f_a[2]; /* adouble for dependent variables */ 68 adouble x_a[2]; /* adouble for independent variables */ 69 70 PetscFunctionBeginUser; 71 PetscCall(VecGetArrayRead(X,&x)); 72 PetscCall(VecGetArray(F,&f)); 73 74 trace_on(1); /* Start of active section */ 75 x_a[0] <<= x[0]; x_a[1] <<= x[1]; /* Mark as independent */ 76 f_a[0] = x_a[1]; 77 f_a[1] = mu*(1.-x_a[0]*x_a[0])*x_a[1]-x_a[0]; 78 f_a[0] >>= f[0]; f_a[1] >>= f[1]; /* Mark as dependent */ 79 trace_off(1); /* End of active section */ 80 81 PetscCall(VecRestoreArrayRead(X,&x)); 82 PetscCall(VecRestoreArray(F,&f)); 83 PetscFunctionReturn(0); 84 } 85 86 /* 87 Compute the Jacobian w.r.t. x using PETSc-ADOL-C driver. 88 */ 89 static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx) 90 { 91 User user=(User)ctx; 92 const PetscScalar *x; 93 94 PetscFunctionBeginUser; 95 PetscCall(VecGetArrayRead(X,&x)); 96 PetscCall(PetscAdolcComputeRHSJacobian(1,A,x,user->adctx)); 97 PetscCall(VecRestoreArrayRead(X,&x)); 98 PetscFunctionReturn(0); 99 } 100 101 /* 102 Monitor timesteps and use interpolation to output at integer multiples of 0.1 103 */ 104 static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 105 { 106 const PetscScalar *x; 107 PetscReal tfinal, dt, tprev; 108 User user = (User)ctx; 109 110 PetscFunctionBeginUser; 111 PetscCall(TSGetTimeStep(ts,&dt)); 112 PetscCall(TSGetMaxTime(ts,&tfinal)); 113 PetscCall(TSGetPrevTime(ts,&tprev)); 114 PetscCall(VecGetArrayRead(X,&x)); 115 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D 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]))); 116 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev)); 117 PetscCall(VecGetArrayRead(X,&x)); 118 PetscFunctionReturn(0); 119 } 120 121 int main(int argc,char **argv) 122 { 123 TS ts = NULL; /* nonlinear solver */ 124 Vec ic,r; 125 PetscBool monitor = PETSC_FALSE; 126 PetscScalar *x_ptr; 127 PetscMPIInt size; 128 struct _n_User user; 129 AdolcCtx *adctx; 130 Tao tao; 131 KSP ksp; 132 PC pc; 133 134 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135 Initialize program 136 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 137 PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 138 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 139 PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 140 141 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142 Set runtime options and create AdolcCtx 143 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144 PetscCall(PetscNew(&adctx)); 145 user.mu = 1.0; 146 user.next_output = 0.0; 147 user.steps = 0; 148 user.ftime = 0.5; 149 adctx->m = 2;adctx->n = 2;adctx->p = 2; 150 user.adctx = adctx; 151 152 PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL)); 153 PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 154 155 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156 Create necessary matrix and vectors, solve same ODE on every process 157 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 158 PetscCall(MatCreate(PETSC_COMM_WORLD,&user.A)); 159 PetscCall(MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2)); 160 PetscCall(MatSetFromOptions(user.A)); 161 PetscCall(MatSetUp(user.A)); 162 PetscCall(MatCreateVecs(user.A,&user.x,NULL)); 163 164 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 165 Set initial conditions 166 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 167 PetscCall(VecGetArray(user.x,&x_ptr)); 168 x_ptr[0] = 2.0; x_ptr[1] = 0.66666654321; 169 PetscCall(VecRestoreArray(user.x,&x_ptr)); 170 171 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 172 Trace just once on each tape and put zeros on Jacobian diagonal 173 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 174 PetscCall(VecDuplicate(user.x,&r)); 175 PetscCall(RHSFunctionActive(ts,0.,user.x,r,&user)); 176 PetscCall(VecSet(r,0)); 177 PetscCall(MatDiagonalSet(user.A,r,INSERT_VALUES)); 178 PetscCall(VecDestroy(&r)); 179 180 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181 Create timestepping solver context 182 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 183 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 184 PetscCall(TSSetType(ts,TSRK)); 185 PetscCall(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&user)); 186 PetscCall(TSSetRHSJacobian(ts,user.A,user.A,RHSJacobian,&user)); 187 PetscCall(TSSetMaxTime(ts,user.ftime)); 188 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 189 if (monitor) { 190 PetscCall(TSMonitorSet(ts,Monitor,&user,NULL)); 191 } 192 193 PetscCall(TSSetTime(ts,0.0)); 194 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)(user.ftime))); 195 196 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 197 Save trajectory of solution so that TSAdjointSolve() may be used 198 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 199 PetscCall(TSSetSaveTrajectory(ts)); 200 201 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 202 Set runtime options 203 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 204 PetscCall(TSSetFromOptions(ts)); 205 206 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 207 Solve nonlinear system 208 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 209 PetscCall(TSSolve(ts,user.x)); 210 PetscCall(TSGetSolveTime(ts,&(user.ftime))); 211 PetscCall(TSGetStepNumber(ts,&user.steps)); 212 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime)); 213 214 PetscCall(VecGetArray(user.x,&x_ptr)); 215 user.x_ob[0] = x_ptr[0]; 216 user.x_ob[1] = x_ptr[1]; 217 PetscCall(VecRestoreArray(user.x,&x_ptr)); 218 219 PetscCall(MatCreateVecs(user.A,&user.lambda[0],NULL)); 220 221 /* Create TAO solver and set desired solution method */ 222 PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao)); 223 PetscCall(TaoSetType(tao,TAOCG)); 224 225 /* Set initial solution guess */ 226 PetscCall(MatCreateVecs(user.A,&ic,NULL)); 227 PetscCall(VecGetArray(ic,&x_ptr)); 228 x_ptr[0] = 2.1; 229 x_ptr[1] = 0.7; 230 PetscCall(VecRestoreArray(ic,&x_ptr)); 231 232 PetscCall(TaoSetSolution(tao,ic)); 233 234 /* Set routine for function and gradient evaluation */ 235 PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user)); 236 237 /* Check for any TAO command line options */ 238 PetscCall(TaoSetFromOptions(tao)); 239 PetscCall(TaoGetKSP(tao,&ksp)); 240 if (ksp) { 241 PetscCall(KSPGetPC(ksp,&pc)); 242 PetscCall(PCSetType(pc,PCNONE)); 243 } 244 245 PetscCall(TaoSetTolerances(tao,1e-10,PETSC_DEFAULT,PETSC_DEFAULT)); 246 247 /* SOLVE THE APPLICATION */ 248 PetscCall(TaoSolve(tao)); 249 250 /* Free TAO data structures */ 251 PetscCall(TaoDestroy(&tao)); 252 253 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 254 Free work space. All PETSc objects should be destroyed when they 255 are no longer needed. 256 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 257 PetscCall(MatDestroy(&user.A)); 258 PetscCall(VecDestroy(&user.x)); 259 PetscCall(VecDestroy(&user.lambda[0])); 260 PetscCall(TSDestroy(&ts)); 261 PetscCall(VecDestroy(&ic)); 262 PetscCall(PetscFree(adctx)); 263 PetscCall(PetscFinalize()); 264 return 0; 265 } 266 267 /* ------------------------------------------------------------------ */ 268 /* 269 FormFunctionGradient - Evaluates the function and corresponding gradient. 270 271 Input Parameters: 272 tao - the Tao context 273 X - the input vector 274 ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 275 276 Output Parameters: 277 f - the newly evaluated function 278 G - the newly evaluated gradient 279 */ 280 PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx) 281 { 282 User user = (User)ctx; 283 TS ts; 284 PetscScalar *x_ptr,*y_ptr; 285 286 PetscFunctionBeginUser; 287 PetscCall(VecCopy(IC,user->x)); 288 289 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 290 Create timestepping solver context 291 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 292 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 293 PetscCall(TSSetType(ts,TSRK)); 294 PetscCall(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,user)); 295 /* Set RHS Jacobian for the adjoint integration */ 296 PetscCall(TSSetRHSJacobian(ts,user->A,user->A,RHSJacobian,user)); 297 298 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 299 Set time 300 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 301 PetscCall(TSSetTime(ts,0.0)); 302 PetscCall(TSSetTimeStep(ts,.001)); 303 PetscCall(TSSetMaxTime(ts,0.5)); 304 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 305 306 PetscCall(TSSetTolerances(ts,1e-7,NULL,1e-7,NULL)); 307 308 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 309 Save trajectory of solution so that TSAdjointSolve() may be used 310 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 311 PetscCall(TSSetSaveTrajectory(ts)); 312 313 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 314 Set runtime options 315 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 316 PetscCall(TSSetFromOptions(ts)); 317 318 PetscCall(TSSolve(ts,user->x)); 319 PetscCall(TSGetSolveTime(ts,&user->ftime)); 320 PetscCall(TSGetStepNumber(ts,&user->steps)); 321 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %.6f, steps %D, ftime %g\n",(double)user->mu,user->steps,(double)user->ftime)); 322 323 PetscCall(VecGetArray(user->x,&x_ptr)); 324 *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]); 325 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))); 326 327 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 328 Adjoint model starts here 329 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 330 /* Redet initial conditions for the adjoint integration */ 331 PetscCall(VecGetArray(user->lambda[0],&y_ptr)); 332 y_ptr[0] = 2.*(x_ptr[0]-user->x_ob[0]); 333 y_ptr[1] = 2.*(x_ptr[1]-user->x_ob[1]); 334 PetscCall(VecRestoreArray(user->lambda[0],&y_ptr)); 335 PetscCall(VecRestoreArray(user->x,&x_ptr)); 336 PetscCall(TSSetCostGradients(ts,1,user->lambda,NULL)); 337 338 PetscCall(TSAdjointSolve(ts)); 339 340 PetscCall(VecCopy(user->lambda[0],G)); 341 342 PetscCall(TSDestroy(&ts)); 343 PetscFunctionReturn(0); 344 } 345 346 /*TEST 347 348 build: 349 requires: double !complex adolc 350 351 test: 352 suffix: 1 353 args: -ts_rhs_jacobian_test_mult_transpose FALSE -tao_max_it 2 -ts_rhs_jacobian_test_mult FALSE 354 output_file: output/ex16opt_ic_1.out 355 356 TEST*/ 357