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