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