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