xref: /petsc/src/ts/tutorials/hybrid/ex1fd.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 static char help[] = "An example of hybrid system using TS event.\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   Reference:
15   I. A. Hiskens, M.A. Pai, Trajectory Sensitivity Analysis of Hybrid Systems, IEEE Transactions on Circuits and Systems, Vol 47, No 2, February 2000
16 */
17 
18 #include <petscts.h>
19 
20 typedef struct {
21   PetscReal lambda1;
22   PetscReal lambda2;
23   PetscInt  mode; /* mode flag*/
24 } AppCtx;
25 
26 PetscErrorCode FWDRun(TS, Vec, void *);
27 
28 PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx) {
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   Vec               *lambda, *mu;
45   PetscScalar       *x, *y;
46   const PetscScalar *u;
47   PetscScalar        tmp[2], A1[2][2], A2[2], denorm1, denorm2;
48   PetscInt           numcost;
49 
50   PetscFunctionBegin;
51   PetscCall(TSGetCostGradients(ts, &numcost, &lambda, &mu));
52   PetscCall(VecGetArrayRead(U, &u));
53 
54   if (actx->mode == 2) {
55     denorm1  = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]);
56     denorm2  = -actx->lambda1 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]);
57     A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm1 + 1.;
58     A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm1;
59     A1[1][0] = 110. * u[1] * 1. / denorm1;
60     A1[1][1] = -110. * u[0] * 1. / denorm1 + 1.;
61 
62     A2[0] = 110. * u[1] * (-u[0]) / denorm2;
63     A2[1] = -110. * u[0] * (-u[0]) / denorm2;
64   } else {
65     denorm2  = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]);
66     A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm2 + 1;
67     A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm2;
68     A1[1][0] = 110. * u[1] * 1. / denorm2;
69     A1[1][1] = -110. * u[0] * 1. / denorm2 + 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   AppCtx *actx = (AppCtx *)ctx;
101 
102   PetscFunctionBegin;
103   if (!forwardsolve) PetscCall(ShiftGradients(ts, U, actx));
104   if (actx->mode == 1) {
105     actx->mode = 2;
106     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Change from mode 1 to 2 at t = %f \n", (double)t));
107   } else if (actx->mode == 2) {
108     actx->mode = 1;
109     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Change from mode 2 to 1 at t = %f \n", (double)t));
110   }
111   PetscFunctionReturn(0);
112 }
113 
114 /*
115      Defines the ODE passed to the ODE solver
116 */
117 static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) {
118   AppCtx            *actx = (AppCtx *)ctx;
119   PetscScalar       *f;
120   const PetscScalar *u, *udot;
121 
122   PetscFunctionBegin;
123   /*  The next three lines allow us to access the entries of the vectors directly */
124   PetscCall(VecGetArrayRead(U, &u));
125   PetscCall(VecGetArrayRead(Udot, &udot));
126   PetscCall(VecGetArray(F, &f));
127 
128   if (actx->mode == 1) {
129     f[0] = udot[0] - u[0] + 100 * u[1];
130     f[1] = udot[1] - 10 * u[0] - u[1];
131   } else if (actx->mode == 2) {
132     f[0] = udot[0] - u[0] - 10 * u[1];
133     f[1] = udot[1] + 100 * u[0] - u[1];
134   }
135 
136   PetscCall(VecRestoreArrayRead(U, &u));
137   PetscCall(VecRestoreArrayRead(Udot, &udot));
138   PetscCall(VecRestoreArray(F, &f));
139   PetscFunctionReturn(0);
140 }
141 
142 /*
143      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
144 */
145 static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) {
146   AppCtx            *actx     = (AppCtx *)ctx;
147   PetscInt           rowcol[] = {0, 1};
148   PetscScalar        J[2][2];
149   const PetscScalar *u, *udot;
150 
151   PetscFunctionBegin;
152   PetscCall(VecGetArrayRead(U, &u));
153   PetscCall(VecGetArrayRead(Udot, &udot));
154 
155   if (actx->mode == 1) {
156     J[0][0] = a - 1;
157     J[0][1] = 100;
158     J[1][0] = -10;
159     J[1][1] = a - 1;
160   } else if (actx->mode == 2) {
161     J[0][0] = a - 1;
162     J[0][1] = -10;
163     J[1][0] = 100;
164     J[1][1] = a - 1;
165   }
166   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
167 
168   PetscCall(VecRestoreArrayRead(U, &u));
169   PetscCall(VecRestoreArrayRead(Udot, &udot));
170 
171   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
172   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
173   if (A != B) {
174     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
175     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
176   }
177   PetscFunctionReturn(0);
178 }
179 
180 int main(int argc, char **argv) {
181   TS           ts; /* ODE integrator */
182   Vec          U;  /* solution will be stored here */
183   Mat          A;  /* Jacobian matrix */
184   Mat          Ap; /* dfdp */
185   PetscMPIInt  size;
186   PetscInt     n = 2;
187   PetscScalar *u;
188   AppCtx       app;
189   PetscInt     direction[1];
190   PetscBool    terminate[1];
191   PetscReal    delta, tmp[2], sensi[2];
192 
193   delta = 1e-8;
194 
195   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196      Initialize program
197      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198   PetscFunctionBeginUser;
199   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
200   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
201   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
202   app.mode    = 1;
203   app.lambda1 = 2.75;
204   app.lambda2 = 0.36;
205   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1 options", "");
206   {
207     PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL));
208     PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL));
209   }
210   PetscOptionsEnd();
211 
212   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
213     Create necessary matrix and vectors
214     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
215   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
216   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
217   PetscCall(MatSetType(A, MATDENSE));
218   PetscCall(MatSetFromOptions(A));
219   PetscCall(MatSetUp(A));
220 
221   PetscCall(MatCreateVecs(A, &U, NULL));
222 
223   PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap));
224   PetscCall(MatSetSizes(Ap, n, 1, PETSC_DETERMINE, PETSC_DETERMINE));
225   PetscCall(MatSetType(Ap, MATDENSE));
226   PetscCall(MatSetFromOptions(Ap));
227   PetscCall(MatSetUp(Ap));
228   PetscCall(MatZeroEntries(Ap)); /* initialize to zeros */
229 
230   PetscCall(VecGetArray(U, &u));
231   u[0] = 0;
232   u[1] = 1;
233   PetscCall(VecRestoreArray(U, &u));
234   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
235      Create timestepping solver context
236      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
237   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
238   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
239   PetscCall(TSSetType(ts, TSCN));
240   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &app));
241   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &app));
242 
243   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244      Set initial conditions
245    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
246   PetscCall(TSSetSolution(ts, U));
247 
248   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
249      Set solver options
250    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
251   PetscCall(TSSetMaxTime(ts, 0.125));
252   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
253   PetscCall(TSSetTimeStep(ts, 1. / 256.));
254   PetscCall(TSSetFromOptions(ts));
255 
256   /* Set directions and terminate flags for the two events */
257   direction[0] = 0;
258   terminate[0] = PETSC_FALSE;
259   PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app));
260 
261   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262      Run timestepping solver
263      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
264   PetscCall(TSSolve(ts, U));
265 
266   PetscCall(VecGetArray(U, &u));
267   tmp[0] = u[0];
268   tmp[1] = u[1];
269 
270   u[0] = 0 + delta;
271   u[1] = 1;
272   PetscCall(VecRestoreArray(U, &u));
273 
274   PetscCall(FWDRun(ts, U, (void *)&app));
275 
276   PetscCall(VecGetArray(U, &u));
277   sensi[0] = (u[0] - tmp[0]) / delta;
278   sensi[1] = (u[1] - tmp[1]) / delta;
279   PetscCall(PetscPrintf(PETSC_COMM_SELF, "d x1(tf) /d x1(t0) = %f d x2(tf) / d x1(t0) = %f \n", (double)sensi[0], (double)sensi[1]));
280   u[0] = 0;
281   u[1] = 1 + delta;
282   PetscCall(VecRestoreArray(U, &u));
283 
284   PetscCall(FWDRun(ts, U, (void *)&app));
285 
286   PetscCall(VecGetArray(U, &u));
287   sensi[0] = (u[0] - tmp[0]) / delta;
288   sensi[1] = (u[1] - tmp[1]) / delta;
289   PetscCall(PetscPrintf(PETSC_COMM_SELF, "d x1(tf) /d x2(t0) = %f d x2(tf) / d x2(t0) = %f \n", (double)sensi[0], (double)sensi[1]));
290   u[0]        = 0;
291   u[1]        = 1;
292   app.lambda1 = app.lambda1 + delta;
293   PetscCall(VecRestoreArray(U, &u));
294 
295   PetscCall(FWDRun(ts, U, (void *)&app));
296 
297   PetscCall(VecGetArray(U, &u));
298   sensi[0] = (u[0] - tmp[0]) / delta;
299   sensi[1] = (u[1] - tmp[1]) / delta;
300   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Final gradients: d x1(tf) /d p = %f d x2(tf) / d p = %f \n", (double)sensi[0], (double)sensi[1]));
301   PetscCall(VecRestoreArray(U, &u));
302 
303   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
304      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
305    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
306   PetscCall(MatDestroy(&A));
307   PetscCall(VecDestroy(&U));
308   PetscCall(TSDestroy(&ts));
309 
310   PetscCall(MatDestroy(&Ap));
311   PetscCall(PetscFinalize());
312   return 0;
313 }
314 
315 PetscErrorCode FWDRun(TS ts, Vec U0, void *ctx0) {
316   Vec     U; /* solution will be stored here */
317   AppCtx *ctx = (AppCtx *)ctx0;
318 
319   PetscFunctionBeginUser;
320   PetscCall(TSGetSolution(ts, &U));
321   PetscCall(VecCopy(U0, U));
322 
323   ctx->mode = 1;
324   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
325      Run timestepping solver
326      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
327   PetscCall(TSSetTime(ts, 0.0));
328 
329   PetscCall(TSSolve(ts, U));
330 
331   PetscFunctionReturn(0);
332 }
333 
334 /*TEST
335 
336    build:
337       requires: !defined(PETSC_USE_CXXCOMPLEX)
338 
339    test:
340       args: -ts_event_tol 1e-9
341       timeoutfactor: 18
342       requires: !single
343 
344 TEST*/
345