xref: /petsc/src/ts/tutorials/hybrid/ex1adj.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
1 static char help[] = "Adjoint sensitivity of a hybrid system with state-dependent switchings.\n";
2 
3 /*
4   The dynamics is described by the ODE
5                   u_t = A_i u
6 
7   where A_1 = [ 1  -100
8                 10  1  ],
9         A_2 = [ 1    10
10                -100  1 ].
11   The index i changes from 1 to 2 when u[1]=2.75u[0] and from 2 to 1 when u[1]=0.36u[0].
12   Initially u=[0 1]^T and i=1.
13 
14   References:
15 + * - H. Zhang, S. Abhyankar, E. Constantinescu, M. Mihai, Discrete Adjoint Sensitivity Analysis of Hybrid Dynamical Systems With Switching, IEEE Transactions on Circuits and Systems I: Regular Papers, 64(5), May 2017
16 - * - I. A. Hiskens, M.A. Pai, Trajectory Sensitivity Analysis of Hybrid Systems, IEEE Transactions on Circuits and Systems, Vol 47, No 2, February 2000
17 */
18 
19 #include <petscts.h>
20 
21 typedef struct {
22   PetscScalar lambda1;
23   PetscScalar lambda2;
24   PetscInt    mode; /* mode flag*/
25 } AppCtx;
26 
27 PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx)
28 {
29   AppCtx            *actx = (AppCtx *)ctx;
30   const PetscScalar *u;
31 
32   PetscFunctionBegin;
33   PetscCall(VecGetArrayRead(U, &u));
34   if (actx->mode == 1) {
35     fvalue[0] = u[1] - actx->lambda1 * u[0];
36   } else if (actx->mode == 2) {
37     fvalue[0] = u[1] - actx->lambda2 * u[0];
38   }
39   PetscCall(VecRestoreArrayRead(U, &u));
40   PetscFunctionReturn(0);
41 }
42 
43 PetscErrorCode ShiftGradients(TS ts, Vec U, AppCtx *actx)
44 {
45   Vec               *lambda, *mu;
46   PetscScalar       *x, *y;
47   const PetscScalar *u;
48   PetscScalar        tmp[2], A1[2][2], A2[2], denorm;
49   PetscInt           numcost;
50 
51   PetscFunctionBegin;
52   PetscCall(TSGetCostGradients(ts, &numcost, &lambda, &mu));
53   PetscCall(VecGetArrayRead(U, &u));
54 
55   if (actx->mode == 2) {
56     denorm   = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]);
57     A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm + 1.;
58     A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm;
59     A1[1][0] = 110. * u[1] * 1. / denorm;
60     A1[1][1] = -110. * u[0] * 1. / denorm + 1.;
61 
62     A2[0] = 110. * u[1] * (-u[0]) / denorm;
63     A2[1] = -110. * u[0] * (-u[0]) / denorm;
64   } else {
65     denorm   = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]);
66     A1[0][0] = 110. * u[1] * (actx->lambda2) / denorm + 1;
67     A1[0][1] = -110. * u[0] * (actx->lambda2) / denorm;
68     A1[1][0] = -110. * u[1] * 1. / denorm;
69     A1[1][1] = 110. * u[0] * 1. / denorm + 1.;
70 
71     A2[0] = 0;
72     A2[1] = 0;
73   }
74 
75   PetscCall(VecRestoreArrayRead(U, &u));
76 
77   PetscCall(VecGetArray(lambda[0], &x));
78   PetscCall(VecGetArray(mu[0], &y));
79   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
80   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
81   y[0]   = y[0] + A2[0] * x[0] + A2[1] * x[1];
82   x[0]   = tmp[0];
83   x[1]   = tmp[1];
84   PetscCall(VecRestoreArray(mu[0], &y));
85   PetscCall(VecRestoreArray(lambda[0], &x));
86 
87   PetscCall(VecGetArray(lambda[1], &x));
88   PetscCall(VecGetArray(mu[1], &y));
89   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
90   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
91   y[0]   = y[0] + A2[0] * x[0] + A2[1] * x[1];
92   x[0]   = tmp[0];
93   x[1]   = tmp[1];
94   PetscCall(VecRestoreArray(mu[1], &y));
95   PetscCall(VecRestoreArray(lambda[1], &x));
96   PetscFunctionReturn(0);
97 }
98 
99 PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx)
100 {
101   AppCtx *actx = (AppCtx *)ctx;
102 
103   PetscFunctionBegin;
104   /* PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); */
105   if (!forwardsolve) PetscCall(ShiftGradients(ts, U, actx));
106   if (actx->mode == 1) {
107     actx->mode = 2;
108     /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t)); */
109   } else if (actx->mode == 2) {
110     actx->mode = 1;
111     /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t)); */
112   }
113   PetscFunctionReturn(0);
114 }
115 
116 /*
117      Defines the ODE passed to the ODE solver
118 */
119 static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx)
120 {
121   AppCtx            *actx = (AppCtx *)ctx;
122   PetscScalar       *f;
123   const PetscScalar *u, *udot;
124 
125   PetscFunctionBegin;
126   /*  The next three lines allow us to access the entries of the vectors directly */
127   PetscCall(VecGetArrayRead(U, &u));
128   PetscCall(VecGetArrayRead(Udot, &udot));
129   PetscCall(VecGetArray(F, &f));
130 
131   if (actx->mode == 1) {
132     f[0] = udot[0] - u[0] + 100 * u[1];
133     f[1] = udot[1] - 10 * u[0] - u[1];
134   } else if (actx->mode == 2) {
135     f[0] = udot[0] - u[0] - 10 * u[1];
136     f[1] = udot[1] + 100 * u[0] - u[1];
137   }
138 
139   PetscCall(VecRestoreArrayRead(U, &u));
140   PetscCall(VecRestoreArrayRead(Udot, &udot));
141   PetscCall(VecRestoreArray(F, &f));
142   PetscFunctionReturn(0);
143 }
144 
145 /*
146      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
147 */
148 static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx)
149 {
150   AppCtx            *actx     = (AppCtx *)ctx;
151   PetscInt           rowcol[] = {0, 1};
152   PetscScalar        J[2][2];
153   const PetscScalar *u, *udot;
154 
155   PetscFunctionBegin;
156   PetscCall(VecGetArrayRead(U, &u));
157   PetscCall(VecGetArrayRead(Udot, &udot));
158 
159   if (actx->mode == 1) {
160     J[0][0] = a - 1;
161     J[0][1] = 100;
162     J[1][0] = -10;
163     J[1][1] = a - 1;
164   } else if (actx->mode == 2) {
165     J[0][0] = a - 1;
166     J[0][1] = -10;
167     J[1][0] = 100;
168     J[1][1] = a - 1;
169   }
170   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
171 
172   PetscCall(VecRestoreArrayRead(U, &u));
173   PetscCall(VecRestoreArrayRead(Udot, &udot));
174 
175   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
176   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
177   if (A != B) {
178     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
179     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
180   }
181   PetscFunctionReturn(0);
182 }
183 
184 /* Matrix JacobianP is constant so that it only needs to be evaluated once */
185 static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx)
186 {
187   PetscFunctionBeginUser;
188   PetscFunctionReturn(0);
189 }
190 
191 int main(int argc, char **argv)
192 {
193   TS           ts; /* ODE integrator */
194   Vec          U;  /* solution will be stored here */
195   Mat          A;  /* Jacobian matrix */
196   Mat          Ap; /* dfdp */
197   PetscMPIInt  size;
198   PetscInt     n = 2;
199   PetscScalar *u, *v;
200   AppCtx       app;
201   PetscInt     direction[1];
202   PetscBool    terminate[1];
203   Vec          lambda[2], mu[2];
204   PetscReal    tend;
205 
206   FILE *f;
207   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208      Initialize program
209      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
210   PetscFunctionBeginUser;
211   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
212   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
213   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
214   app.mode    = 1;
215   app.lambda1 = 2.75;
216   app.lambda2 = 0.36;
217   tend        = 0.125;
218   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1adj options", "");
219   {
220     PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL));
221     PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL));
222     PetscCall(PetscOptionsReal("-tend", "", "", tend, &tend, NULL));
223   }
224   PetscOptionsEnd();
225 
226   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227     Create necessary matrix and vectors
228     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
230   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
231   PetscCall(MatSetType(A, MATDENSE));
232   PetscCall(MatSetFromOptions(A));
233   PetscCall(MatSetUp(A));
234 
235   PetscCall(MatCreateVecs(A, &U, NULL));
236 
237   PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap));
238   PetscCall(MatSetSizes(Ap, n, 1, PETSC_DETERMINE, PETSC_DETERMINE));
239   PetscCall(MatSetType(Ap, MATDENSE));
240   PetscCall(MatSetFromOptions(Ap));
241   PetscCall(MatSetUp(Ap));
242   PetscCall(MatZeroEntries(Ap)); /* initialize to zeros */
243 
244   PetscCall(VecGetArray(U, &u));
245   u[0] = 0;
246   u[1] = 1;
247   PetscCall(VecRestoreArray(U, &u));
248   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
249      Create timestepping solver context
250      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
251   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
252   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
253   PetscCall(TSSetType(ts, TSCN));
254   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &app));
255   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &app));
256   PetscCall(TSSetRHSJacobianP(ts, Ap, RHSJacobianP, &app));
257 
258   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
259      Set initial conditions
260    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
261   PetscCall(TSSetSolution(ts, U));
262 
263   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
264     Save trajectory of solution so that TSAdjointSolve() may be used
265    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
266   PetscCall(TSSetSaveTrajectory(ts));
267 
268   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269      Set solver options
270    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
271   PetscCall(TSSetMaxTime(ts, tend));
272   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
273   PetscCall(TSSetTimeStep(ts, 1. / 256.));
274   PetscCall(TSSetFromOptions(ts));
275 
276   /* Set directions and terminate flags for the two events */
277   direction[0] = 0;
278   terminate[0] = PETSC_FALSE;
279   PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app));
280 
281   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
282      Run timestepping solver
283      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
284   PetscCall(TSSolve(ts, U));
285 
286   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
287      Adjoint model starts here
288      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
289   PetscCall(MatCreateVecs(A, &lambda[0], NULL));
290   PetscCall(MatCreateVecs(A, &lambda[1], NULL));
291   /*   Set initial conditions for the adjoint integration */
292   PetscCall(VecZeroEntries(lambda[0]));
293   PetscCall(VecZeroEntries(lambda[1]));
294   PetscCall(VecGetArray(lambda[0], &u));
295   u[0] = 1.;
296   PetscCall(VecRestoreArray(lambda[0], &u));
297   PetscCall(VecGetArray(lambda[1], &u));
298   u[1] = 1.;
299   PetscCall(VecRestoreArray(lambda[1], &u));
300 
301   PetscCall(MatCreateVecs(Ap, &mu[0], NULL));
302   PetscCall(MatCreateVecs(Ap, &mu[1], NULL));
303   PetscCall(VecZeroEntries(mu[0]));
304   PetscCall(VecZeroEntries(mu[1]));
305   PetscCall(TSSetCostGradients(ts, 2, lambda, mu));
306 
307   PetscCall(TSAdjointSolve(ts));
308 
309   /*
310   PetscCall(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD));
311   PetscCall(VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD));
312   PetscCall(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD));
313   PetscCall(VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD));
314   */
315   PetscCall(VecGetArray(mu[0], &u));
316   PetscCall(VecGetArray(mu[1], &v));
317   f = fopen("adj_mu.out", "a");
318   PetscCall(PetscFPrintf(PETSC_COMM_WORLD, f, "%20.15lf %20.15lf %20.15lf\n", (double)tend, (double)PetscRealPart(u[0]), (double)PetscRealPart(v[0])));
319   PetscCall(VecRestoreArray(mu[0], &u));
320   PetscCall(VecRestoreArray(mu[1], &v));
321   fclose(f);
322   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
323      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
324    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
325   PetscCall(MatDestroy(&A));
326   PetscCall(VecDestroy(&U));
327   PetscCall(TSDestroy(&ts));
328 
329   PetscCall(MatDestroy(&Ap));
330   PetscCall(VecDestroy(&lambda[0]));
331   PetscCall(VecDestroy(&lambda[1]));
332   PetscCall(VecDestroy(&mu[0]));
333   PetscCall(VecDestroy(&mu[1]));
334   PetscCall(PetscFinalize());
335   return 0;
336 }
337 
338 /*TEST
339 
340    build:
341       requires: !complex
342 
343    test:
344       args: -ts_monitor -ts_adjoint_monitor
345 
346 TEST*/
347