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