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