xref: /petsc/src/ts/tutorials/autodiff/ex16adj.cxx (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
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   PetscCall(VecGetArrayRead(X,&x));
49   PetscCall(VecGetArray(F,&f));
50   f[0] = x[1];
51   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
52   PetscCall(VecRestoreArrayRead(X,&x));
53   PetscCall(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   PetscCall(VecGetArrayRead(X,&x));
72   PetscCall(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   PetscCall(VecRestoreArrayRead(X,&x));
84   PetscCall(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   PetscCall(VecGetArrayRead(X,&x));
103   PetscCall(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   PetscCall(VecRestoreArrayRead(X,&x));
115   PetscCall(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   PetscCall(VecGetArrayRead(X,&x));
129   PetscCall(PetscAdolcComputeRHSJacobian(1,A,x,user->adctx));
130   PetscCall(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   PetscCall(VecGetArrayRead(X,&x));
144   PetscCall(PetscAdolcComputeRHSJacobianP(3,A,x,&user->mu,user->adctx));
145   PetscCall(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   PetscCall(TSGetTimeStep(ts,&dt));
160   PetscCall(TSGetMaxTime(ts,&tfinal));
161   PetscCall(TSGetPrevTime(ts,&tprev));
162   PetscCall(VecGetArrayRead(X,&x));
163   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])));
164   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev));
165   PetscCall(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   Vec            lambda[2],mu[2],r;
183 
184   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185      Initialize program
186      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
188   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
189   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
190 
191   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192     Set runtime options and create AdolcCtx
193     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
194   PetscCall(PetscNew(&adctx));
195   user.mu          = 1;
196   user.next_output = 0.0;
197   adctx->m = 2;adctx->n = 2;adctx->p = 2;adctx->num_params = 1;
198   user.adctx = adctx;
199 
200   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
201   PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
202 
203   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204     Create necessary matrix and vectors, solve same ODE on every process
205     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
207   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2));
208   PetscCall(MatSetFromOptions(A));
209   PetscCall(MatSetUp(A));
210   PetscCall(MatCreateVecs(A,&x,NULL));
211 
212   PetscCall(MatCreate(PETSC_COMM_WORLD,&Jacp));
213   PetscCall(MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1));
214   PetscCall(MatSetFromOptions(Jacp));
215   PetscCall(MatSetUp(Jacp));
216 
217   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218      Create timestepping solver context
219      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
220   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
221   PetscCall(TSSetType(ts,TSRK));
222   PetscCall(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&user));
223 
224   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225      Set initial conditions
226    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227   PetscCall(VecGetArray(x,&x_ptr));
228   x_ptr[0] = 2;   x_ptr[1] = 0.66666654321;
229   PetscCall(VecRestoreArray(x,&x_ptr));
230 
231   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232      Trace just once on each tape and put zeros on Jacobian diagonal
233      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
234   PetscCall(VecDuplicate(x,&r));
235   PetscCall(RHSFunctionActive(ts,0.,x,r,&user));
236   PetscCall(RHSFunctionActiveP(ts,0.,x,r,&user));
237   PetscCall(VecSet(r,0));
238   PetscCall(MatDiagonalSet(A,r,INSERT_VALUES));
239   PetscCall(VecDestroy(&r));
240 
241   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
242      Set RHS Jacobian for the adjoint integration
243      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
244   PetscCall(TSSetRHSJacobian(ts,A,A,RHSJacobian,&user));
245   PetscCall(TSSetMaxTime(ts,ftime));
246   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
247   if (monitor) {
248     PetscCall(TSMonitorSet(ts,Monitor,&user,NULL));
249   }
250   PetscCall(TSSetTimeStep(ts,.001));
251 
252   /*
253     Have the TS save its trajectory so that TSAdjointSolve() may be used
254   */
255   PetscCall(TSSetSaveTrajectory(ts));
256 
257   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258      Set runtime options
259    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
260   PetscCall(TSSetFromOptions(ts));
261 
262   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
263      Solve nonlinear system
264      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
265   PetscCall(TSSolve(ts,x));
266   PetscCall(TSGetSolveTime(ts,&ftime));
267   PetscCall(TSGetStepNumber(ts,&steps));
268   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime));
269   PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
270 
271   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272      Start the Adjoint model
273      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
274   PetscCall(MatCreateVecs(A,&lambda[0],NULL));
275   PetscCall(MatCreateVecs(A,&lambda[1],NULL));
276   /*   Reset initial conditions for the adjoint integration */
277   PetscCall(VecGetArray(lambda[0],&x_ptr));
278   x_ptr[0] = 1.0;   x_ptr[1] = 0.0;
279   PetscCall(VecRestoreArray(lambda[0],&x_ptr));
280   PetscCall(VecGetArray(lambda[1],&x_ptr));
281   x_ptr[0] = 0.0;   x_ptr[1] = 1.0;
282   PetscCall(VecRestoreArray(lambda[1],&x_ptr));
283 
284   PetscCall(MatCreateVecs(Jacp,&mu[0],NULL));
285   PetscCall(MatCreateVecs(Jacp,&mu[1],NULL));
286   PetscCall(VecGetArray(mu[0],&x_ptr));
287   x_ptr[0] = 0.0;
288   PetscCall(VecRestoreArray(mu[0],&x_ptr));
289   PetscCall(VecGetArray(mu[1],&x_ptr));
290   x_ptr[0] = 0.0;
291   PetscCall(VecRestoreArray(mu[1],&x_ptr));
292   PetscCall(TSSetCostGradients(ts,2,lambda,mu));
293 
294   /*   Set RHS JacobianP */
295   PetscCall(TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user));
296 
297   PetscCall(TSAdjointSolve(ts));
298 
299   PetscCall(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD));
300   PetscCall(VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD));
301   PetscCall(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD));
302   PetscCall(VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD));
303 
304   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
305      Free work space.  All PETSc objects should be destroyed when they
306      are no longer needed.
307    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
308   PetscCall(MatDestroy(&A));
309   PetscCall(MatDestroy(&Jacp));
310   PetscCall(VecDestroy(&x));
311   PetscCall(VecDestroy(&lambda[0]));
312   PetscCall(VecDestroy(&lambda[1]));
313   PetscCall(VecDestroy(&mu[0]));
314   PetscCall(VecDestroy(&mu[1]));
315   PetscCall(TSDestroy(&ts));
316   PetscCall(PetscFree(adctx));
317   PetscCall(PetscFinalize());
318   return 0;
319 }
320 
321 /*TEST
322 
323   build:
324     requires: double !complex adolc
325 
326   test:
327     suffix: 1
328     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor
329     output_file: output/ex16adj_1.out
330 
331   test:
332     suffix: 2
333     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -mu 5
334     output_file: output/ex16adj_2.out
335 
336   test:
337     suffix: 3
338     args: -ts_max_steps 10 -monitor
339     output_file: output/ex16adj_3.out
340 
341   test:
342     suffix: 4
343     args: -ts_max_steps 10 -monitor -mu 5
344     output_file: output/ex16adj_4.out
345 
346 TEST*/
347