xref: /petsc/src/ts/tutorials/hybrid/ex1fwd.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 static char help[] = "Trajectory 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   PetscReal   print_time;
26 } AppCtx;
27 
28 PetscErrorCode MyMonitor(TS ts, PetscInt stepnum, PetscReal time, Vec U, void *ctx) {
29   AppCtx      *actx = (AppCtx *)ctx;
30   Mat          sp;
31   PetscScalar *u;
32   PetscInt     nump;
33   FILE        *f;
34 
35   PetscFunctionBegin;
36   if (time >= actx->print_time) {
37     actx->print_time += 1. / 256.;
38     PetscCall(TSForwardGetSensitivities(ts, &nump, &sp));
39     PetscCall(MatDenseGetColumn(sp, 2, &u));
40     f = fopen("fwd_sp.out", "a");
41     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, f, "%20.15lf %20.15lf %20.15lf\n", (double)time, (double)PetscRealPart(u[0]), (double)PetscRealPart(u[1])));
42     PetscCall(MatDenseRestoreColumn(sp, &u));
43     fclose(f);
44   }
45   PetscFunctionReturn(0);
46 }
47 
48 PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx) {
49   AppCtx            *actx = (AppCtx *)ctx;
50   const PetscScalar *u;
51 
52   PetscFunctionBegin;
53   PetscCall(VecGetArrayRead(U, &u));
54   if (actx->mode == 1) {
55     fvalue[0] = u[1] - actx->lambda1 * u[0];
56   } else if (actx->mode == 2) {
57     fvalue[0] = u[1] - actx->lambda2 * u[0];
58   }
59   PetscCall(VecRestoreArrayRead(U, &u));
60   PetscFunctionReturn(0);
61 }
62 
63 PetscErrorCode ShiftGradients(TS ts, Vec U, AppCtx *actx) {
64   Mat          sp;
65   PetscScalar *x;
66   PetscScalar *u;
67   PetscScalar  tmp[2], A1[2][2], A2[2], denorm;
68   PetscInt     nump;
69 
70   PetscFunctionBegin;
71   PetscCall(TSForwardGetSensitivities(ts, &nump, &sp));
72   PetscCall(VecGetArray(U, &u));
73 
74   if (actx->mode == 1) {
75     denorm   = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]);
76     A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm + 1.;
77     A1[1][0] = -110. * u[0] * (-actx->lambda1) / denorm;
78     A1[0][1] = 110. * u[1] * 1. / denorm;
79     A1[1][1] = -110. * u[0] * 1. / denorm + 1.;
80 
81     A2[0] = 110. * u[1] * (-u[0]) / denorm;
82     A2[1] = -110. * u[0] * (-u[0]) / denorm;
83   } else {
84     denorm   = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]);
85     A1[0][0] = 110. * u[1] * (actx->lambda2) / denorm + 1;
86     A1[1][0] = -110. * u[0] * (actx->lambda2) / denorm;
87     A1[0][1] = -110. * u[1] * 1. / denorm;
88     A1[1][1] = 110. * u[0] * 1. / denorm + 1.;
89 
90     A2[0] = 0;
91     A2[1] = 0;
92   }
93 
94   PetscCall(VecRestoreArray(U, &u));
95 
96   PetscCall(MatDenseGetColumn(sp, 0, &x));
97   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
98   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
99   x[0]   = tmp[0];
100   x[1]   = tmp[1];
101   PetscCall(MatDenseRestoreColumn(sp, &x));
102 
103   PetscCall(MatDenseGetColumn(sp, 1, &x));
104   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
105   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
106   x[0]   = tmp[0];
107   x[1]   = tmp[1];
108   PetscCall(MatDenseRestoreColumn(sp, &x));
109 
110   PetscCall(MatDenseGetColumn(sp, 2, &x));
111   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
112   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
113   x[0]   = tmp[0] + A2[0];
114   x[1]   = tmp[1] + A2[1];
115   PetscCall(MatDenseRestoreColumn(sp, &x));
116 
117   PetscFunctionReturn(0);
118 }
119 
120 PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx) {
121   AppCtx *actx = (AppCtx *)ctx;
122 
123   PetscFunctionBegin;
124   /* PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); */
125   PetscCall(ShiftGradients(ts, U, actx));
126   if (actx->mode == 1) {
127     actx->mode = 2;
128     /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t)); */
129   } else if (actx->mode == 2) {
130     actx->mode = 1;
131     /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t)); */
132   }
133   PetscFunctionReturn(0);
134 }
135 
136 /*
137      Defines the ODE passed to the ODE solver
138 */
139 static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) {
140   AppCtx            *actx = (AppCtx *)ctx;
141   PetscScalar       *f;
142   const PetscScalar *u, *udot;
143 
144   PetscFunctionBegin;
145   /*  The next three lines allow us to access the entries of the vectors directly */
146   PetscCall(VecGetArrayRead(U, &u));
147   PetscCall(VecGetArrayRead(Udot, &udot));
148   PetscCall(VecGetArray(F, &f));
149 
150   if (actx->mode == 1) {
151     f[0] = udot[0] - u[0] + 100 * u[1];
152     f[1] = udot[1] - 10 * u[0] - u[1];
153   } else if (actx->mode == 2) {
154     f[0] = udot[0] - u[0] - 10 * u[1];
155     f[1] = udot[1] + 100 * u[0] - u[1];
156   }
157 
158   PetscCall(VecRestoreArrayRead(U, &u));
159   PetscCall(VecRestoreArrayRead(Udot, &udot));
160   PetscCall(VecRestoreArray(F, &f));
161   PetscFunctionReturn(0);
162 }
163 
164 /*
165      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
166 */
167 static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) {
168   AppCtx            *actx     = (AppCtx *)ctx;
169   PetscInt           rowcol[] = {0, 1};
170   PetscScalar        J[2][2];
171   const PetscScalar *u, *udot;
172 
173   PetscFunctionBegin;
174   PetscCall(VecGetArrayRead(U, &u));
175   PetscCall(VecGetArrayRead(Udot, &udot));
176 
177   if (actx->mode == 1) {
178     J[0][0] = a - 1;
179     J[0][1] = 100;
180     J[1][0] = -10;
181     J[1][1] = a - 1;
182   } else if (actx->mode == 2) {
183     J[0][0] = a - 1;
184     J[0][1] = -10;
185     J[1][0] = 100;
186     J[1][1] = a - 1;
187   }
188   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
189 
190   PetscCall(VecRestoreArrayRead(U, &u));
191   PetscCall(VecRestoreArrayRead(Udot, &udot));
192 
193   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
194   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
195   if (A != B) {
196     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
197     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
198   }
199   PetscFunctionReturn(0);
200 }
201 
202 /* Matrix JacobianP is constant so that it only needs to be evaluated once */
203 static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat Ap, void *ctx) {
204   PetscFunctionBeginUser;
205   PetscFunctionReturn(0);
206 }
207 
208 int main(int argc, char **argv) {
209   TS           ts; /* ODE integrator */
210   Vec          U;  /* solution will be stored here */
211   Mat          A;  /* Jacobian matrix */
212   Mat          Ap; /* Jacobian dfdp */
213   PetscMPIInt  size;
214   PetscInt     n = 2;
215   PetscScalar *u;
216   AppCtx       app;
217   PetscInt     direction[1];
218   PetscBool    terminate[1];
219   Mat          sp;
220   PetscReal    tend;
221   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222      Initialize program
223      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224   PetscFunctionBeginUser;
225   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
226   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
227   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
228   app.mode       = 1;
229   app.lambda1    = 2.75;
230   app.lambda2    = 0.36;
231   app.print_time = 1. / 256.;
232   tend           = 0.125;
233   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1fwd options", "");
234   {
235     PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL));
236     PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL));
237     PetscCall(PetscOptionsReal("-tend", "", "", tend, &tend, NULL));
238   }
239   PetscOptionsEnd();
240 
241   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
242     Create necessary matrix and vectors
243     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
244   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
245   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
246   PetscCall(MatSetType(A, MATDENSE));
247   PetscCall(MatSetFromOptions(A));
248   PetscCall(MatSetUp(A));
249 
250   PetscCall(MatCreateVecs(A, &U, NULL));
251 
252   PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap));
253   PetscCall(MatSetSizes(Ap, n, 3, PETSC_DETERMINE, PETSC_DETERMINE));
254   PetscCall(MatSetType(Ap, MATDENSE));
255   PetscCall(MatSetFromOptions(Ap));
256   PetscCall(MatSetUp(Ap));
257   PetscCall(MatZeroEntries(Ap));
258 
259   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, n, 3, NULL, &sp));
260   PetscCall(MatZeroEntries(sp));
261   PetscCall(MatShift(sp, 1.0));
262 
263   PetscCall(VecGetArray(U, &u));
264   u[0] = 0;
265   u[1] = 1;
266   PetscCall(VecRestoreArray(U, &u));
267 
268   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269      Create timestepping solver context
270      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
271   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
272   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
273   PetscCall(TSSetType(ts, TSCN));
274   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &app));
275   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &app));
276 
277   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
278      Set initial conditions
279    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
280   PetscCall(TSSetSolution(ts, U));
281 
282   PetscCall(TSForwardSetSensitivities(ts, 3, sp));
283   /*   Set RHS JacobianP */
284   PetscCall(TSSetRHSJacobianP(ts, Ap, RHSJacobianP, &app));
285 
286   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
287      Set solver options
288    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
289   PetscCall(TSSetMaxTime(ts, tend));
290   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
291   PetscCall(TSSetTimeStep(ts, 1. / 256.));
292   PetscCall(TSMonitorSet(ts, MyMonitor, &app, PETSC_NULL));
293   PetscCall(TSSetFromOptions(ts));
294 
295   /* Set directions and terminate flags for the two events */
296   direction[0] = 0;
297   terminate[0] = PETSC_FALSE;
298   PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app));
299   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
300      Run timestepping solver
301      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
302   PetscCall(TSSolve(ts, U));
303 
304   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
305      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
306    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
307   PetscCall(MatDestroy(&A));
308   PetscCall(VecDestroy(&U));
309   PetscCall(TSDestroy(&ts));
310 
311   PetscCall(MatDestroy(&Ap));
312   PetscCall(MatDestroy(&sp));
313   PetscCall(PetscFinalize());
314   return (0);
315 }
316 
317 /*TEST
318 
319    build:
320       requires: !complex
321 
322    test:
323       args: -ts_monitor
324 
325 TEST*/
326