Actual source code: ex20adj.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: #define c11 1.0
  2: #define c12 0
  3: #define c21 2.0
  4: #define c22 1.0
  5: static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n";

  7: /*
  8:    Concepts: TS^time-dependent nonlinear problems
  9:    Concepts: TS^van der Pol equation DAE equivalent
 10:    Concepts: TS^adjoint sensitivity analysis
 11:    Processors: 1
 12: */
 13: /* ------------------------------------------------------------------------

 15:    This program solves the van der Pol DAE ODE equivalent
 16:       [ u_1' ] = [      u_2              ]  (2)
 17:       [ u_2' ]   [ mu[(1-u_1^2)u_2-u_1]  ]
 18:    on the domain 0 <= x <= 1, with the boundary conditions
 19:        u_1(0) = 2, u_2(0) = -6.666665432100101e-01,
 20:    and
 21:        mu = 10^6,
 22:    and computes the sensitivities of the final solution w.r.t. initial conditions and parameter \mu with the implicit theta method and its discrete adjoint.

 24:    Notes:
 25:    This code demonstrates the TSAdjoint interface to a DAE system.

 27:    The user provides the implicit right-hand-side function
 28:    [ G(u',u,t) ] = [u' - f(u,t)] = [ u_1'] - [        u_2           ]
 29:                                    [ u_2']   [ mu[(1-u_1^2)u_2-u_1] ]

 31:    and the Jacobian of G (from the PETSc user manual)

 33:               dG   dG
 34:    J(G) = a * -- + --
 35:               du'  du

 37:    and the JacobianP of the explicit right-hand side of (2) f(u,t) ( which is equivalent to -G(0,u,t) ).
 38:    df   [       0         ]
 39:    -- = [                 ]
 40:    dp   [ (1 - u_1^2) u_2 ].

 42:    See ex20.c for more details on the Jacobian.

 44:    Many DAEs can be represented in a general form M u_t = f(u,t).
 45:    Thus both sides of (1) are multiplied by an artificial matrix
 46:    M = [ c11 c12 ]
 47:        [ c21 c22 ]
 48:    to turn (1) into the general form. This operation does not change the solution and it is intended for illustration only.

 50:   ------------------------------------------------------------------------- */
 51:  #include <petscts.h>
 52:  #include <petsctao.h>

 54: typedef struct _n_User *User;
 55: struct _n_User {
 56:   PetscReal mu;
 57:   PetscReal next_output;

 59:   /* Sensitivity analysis support */
 60:   PetscInt  steps;
 61:   PetscReal ftime;
 62:   Mat       A;                       /* Jacobian matrix */
 63:   Mat       Jacp;                    /* JacobianP matrix */
 64:   Vec       x,lambda[2],mup[2];  /* adjoint variables */
 65: };

 67: /*
 68: *  User-defined routines
 69: */
 70: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 71: {
 72:   PetscErrorCode    ierr;
 73:   User              user = (User)ctx;
 74:   const PetscScalar *x,*xdot;
 75:   PetscScalar       *f;

 78:   VecGetArrayRead(X,&x);
 79:   VecGetArrayRead(Xdot,&xdot);
 80:   VecGetArray(F,&f);
 81:   f[0] = xdot[0] - x[1];
 82:   f[1] = c21*(xdot[0]-x[1]) + xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]) ;
 83:   VecRestoreArrayRead(X,&x);
 84:   VecRestoreArrayRead(Xdot,&xdot);
 85:   VecRestoreArray(F,&f);
 86:   return(0);
 87: }

 89: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
 90: {
 91:   PetscErrorCode    ierr;
 92:   User              user     = (User)ctx;
 93:   PetscInt          rowcol[] = {0,1};
 94:   PetscScalar       J[2][2];
 95:   const PetscScalar *x;

 98:   VecGetArrayRead(X,&x);

100:   J[0][0] = a;     J[0][1] =  -1.0;
101:   J[1][0] = c21*a + user->mu*(1.0 + 2.0*x[0]*x[1]);   J[1][1] = -c21 + a - user->mu*(1.0-x[0]*x[0]);

103:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
104:   VecRestoreArrayRead(X,&x);

106:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
107:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
108:   if (A != B) {
109:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
110:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
111:   }
112:   return(0);
113: }

115: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
116: {
117:   PetscErrorCode    ierr;
118:   PetscInt          row[] = {0,1},col[]={0};
119:   PetscScalar       J[2][1];
120:   const PetscScalar *x;

123:   VecGetArrayRead(X,&x);

125:   J[0][0] = 0;
126:   J[1][0] = (1.-x[0]*x[0])*x[1]-x[0];
127:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);

129:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
130:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
131:   return(0);
132: }

134: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
135: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
136: {
137:   PetscErrorCode    ierr;
138:   const PetscScalar *x;
139:   PetscReal         tfinal, dt;
140:   User              user = (User)ctx;
141:   Vec               interpolatedX;

144:   TSGetTimeStep(ts,&dt);
145:   TSGetMaxTime(ts,&tfinal);

147:   while (user->next_output <= t && user->next_output <= tfinal) {
148:     VecDuplicate(X,&interpolatedX);
149:     TSInterpolate(ts,user->next_output,interpolatedX);
150:     VecGetArrayRead(interpolatedX,&x);
151:     PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
152:                        user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(x[0]),
153:                        (double)PetscRealPart(x[1]));
154:     VecRestoreArrayRead(interpolatedX,&x);
155:     VecDestroy(&interpolatedX);
156:     user->next_output += 0.1;
157:   }
158:   return(0);
159: }

161: int main(int argc,char **argv)
162: {
163:   TS             ts;            /* nonlinear solver */
164:   PetscBool      monitor = PETSC_FALSE;
165:   PetscScalar    *x_ptr,*y_ptr;
166:   PetscMPIInt    size;
167:   struct _n_User user;

170:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171:      Initialize program
172:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173:   PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
174:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
175:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

177:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178:     Set runtime options
179:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180:   user.next_output = 0.0;
181:   user.mu          = 1.0e6;
182:   user.steps       = 0;
183:   user.ftime       = 0.5;
184:   PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
185:   PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);

187:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188:     Create necessary matrix and vectors, solve same ODE on every process
189:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
190:   MatCreate(PETSC_COMM_WORLD,&user.A);
191:   MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);
192:   MatSetFromOptions(user.A);
193:   MatSetUp(user.A);
194:   MatCreateVecs(user.A,&user.x,NULL);

196:   MatCreate(PETSC_COMM_WORLD,&user.Jacp);
197:   MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
198:   MatSetFromOptions(user.Jacp);
199:   MatSetUp(user.Jacp);

201:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202:      Create timestepping solver context
203:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
204:   TSCreate(PETSC_COMM_WORLD,&ts);
205:   TSSetType(ts,TSCN);
206:   TSSetIFunction(ts,NULL,IFunction,&user);
207:   TSSetIJacobian(ts,user.A,user.A,IJacobian,&user);
208:   TSSetMaxTime(ts,user.ftime);
209:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
210:   if (monitor) {
211:     TSMonitorSet(ts,Monitor,&user,NULL);
212:   }

214:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215:      Set initial conditions
216:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217:   VecGetArray(user.x,&x_ptr);
218:   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
219:   VecRestoreArray(user.x,&x_ptr);
220:   TSSetTimeStep(ts,.0001);

222:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223:     Save trajectory of solution so that TSAdjointSolve() may be used
224:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
225:   TSSetSaveTrajectory(ts);

227:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228:      Set runtime options
229:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
230:   TSSetFromOptions(ts);

232:   TSSolve(ts,user.x);
233:   TSGetSolveTime(ts,&user.ftime);
234:   TSGetStepNumber(ts,&user.steps);

236:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237:      Adjoint model starts here
238:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
239:   MatCreateVecs(user.A,&user.lambda[0],NULL);
240:   /*   Set initial conditions for the adjoint integration */
241:   VecGetArray(user.lambda[0],&y_ptr);
242:   y_ptr[0] = 1.0; y_ptr[1] = 0.0;
243:   VecRestoreArray(user.lambda[0],&y_ptr);
244:   MatCreateVecs(user.A,&user.lambda[1],NULL);
245:   VecGetArray(user.lambda[1],&y_ptr);
246:   y_ptr[0] = 0.0; y_ptr[1] = 1.0;
247:   VecRestoreArray(user.lambda[1],&y_ptr);

249:   MatCreateVecs(user.Jacp,&user.mup[0],NULL);
250:   VecGetArray(user.mup[0],&x_ptr);
251:   x_ptr[0] = 0.0;
252:   VecRestoreArray(user.mup[0],&x_ptr);
253:   MatCreateVecs(user.Jacp,&user.mup[1],NULL);
254:   VecGetArray(user.mup[1],&x_ptr);
255:   x_ptr[0] = 0.0;
256:   VecRestoreArray(user.mup[1],&x_ptr);

258:   TSSetCostGradients(ts,2,user.lambda,user.mup);

260:   /*   Set RHS JacobianP */
261:   TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP,&user);

263:   TSAdjointSolve(ts);

265:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[y(tf)]/d[y0]  d[y(tf)]/d[z0]\n");
266:   VecView(user.lambda[0],PETSC_VIEWER_STDOUT_WORLD);
267:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[z(tf)]/d[y0]  d[z(tf)]/d[z0]\n");
268:   VecView(user.lambda[1],PETSC_VIEWER_STDOUT_WORLD);
269:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameters: d[y(tf)]/d[mu]\n");
270:   VecView(user.mup[0],PETSC_VIEWER_STDOUT_WORLD);
271:   PetscPrintf(PETSC_COMM_WORLD,"\n sensivitity wrt parameters: d[z(tf)]/d[mu]\n");
272:   VecView(user.mup[1],PETSC_VIEWER_STDOUT_WORLD);

274:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
275:      Free work space.  All PETSc objects should be destroyed when they
276:      are no longer needed.
277:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
278:   MatDestroy(&user.A);
279:   MatDestroy(&user.Jacp);
280:   VecDestroy(&user.x);
281:   VecDestroy(&user.lambda[0]);
282:   VecDestroy(&user.lambda[1]);
283:   VecDestroy(&user.mup[0]);
284:   VecDestroy(&user.mup[1]);
285:   TSDestroy(&ts);

287:   PetscFinalize();
288:   return(ierr);
289: }

291: /*TEST

293:     test:
294:       requires: revolve
295:       args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -viewer_binary_skip_info -ts_dt 0.001 -mu 100000 -ts_trajectory_dirname ex20adj1dir

297:     test:
298:       suffix: 2
299:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only

301:     test:
302:       suffix: 3
303:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 0
304:       output_file: output/ex20adj_2.out

306:     test:
307:       suffix: 4
308:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack
309:       output_file: output/ex20adj_2.out

311:     test:
312:       suffix: 5
313:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack
314:       output_file: output/ex20adj_2.out

316:     test:
317:       suffix: 6
318:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack 0
319:       output_file: output/ex20adj_2.out

321:     test:
322:       suffix: 7
323:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0
324:       output_file: output/ex20adj_2.out

326:     test:
327:       suffix: 8
328:       requires: revolve
329:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_solution_only -ts_trajectory_monitor
330:       output_file: output/ex20adj_3.out

332:     test:
333:       suffix: 9
334:       requires: revolve
335:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_solution_only 0 -ts_trajectory_monitor
336:       output_file: output/ex20adj_4.out

338:     test:
339:       requires: revolve
340:       suffix: 10
341:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_revolve_online -ts_trajectory_solution_only
342:       output_file: output/ex20adj_2.out

344:     test:
345:       requires: revolve
346:       suffix: 11
347:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_revolve_online -ts_trajectory_solution_only 0
348:       output_file: output/ex20adj_2.out

350:     test:
351:       suffix: 12
352:       requires: revolve
353:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only
354:       output_file: output/ex20adj_2.out

356:     test:
357:       suffix: 13
358:       requires: revolve
359:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only 0
360:       output_file: output/ex20adj_2.out

362:     test:
363:       suffix: 14
364:       requires: revolve
365:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack
366:       output_file: output/ex20adj_2.out

368:     test:
369:       suffix: 15
370:       requires: revolve
371:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack 0
372:       output_file: output/ex20adj_2.out

374:     test:
375:       suffix: 16
376:       requires: revolve
377:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack
378:       output_file: output/ex20adj_2.out

380:     test:
381:       suffix: 17
382:       requires: revolve
383:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0
384:       output_file: output/ex20adj_2.out

386:     test:
387:       suffix: 18
388:       requires: revolve
389:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack
390:       output_file: output/ex20adj_2.out

392:     test:
393:       suffix: 19
394:       requires: revolve
395:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack
396:       output_file: output/ex20adj_2.out

398:     test:
399:       suffix: 20
400:       requires: revolve
401:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only 0
402:       output_file: output/ex20adj_2.out

404:     test:
405:       suffix: 21
406:       requires: revolve
407:       args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0
408:       output_file: output/ex20adj_2.out

410: TEST*/