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