xref: /petsc/src/ts/tutorials/hybrid/ex1.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
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   PetscScalar lambda1;
22   PetscScalar lambda2;
23   PetscInt    mode;  /* mode flag*/
24 } AppCtx;
25 
26 PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx)
27 {
28   AppCtx         *actx=(AppCtx*)ctx;
29   const PetscScalar *u;
30 
31   PetscFunctionBegin;
32   PetscCall(VecGetArrayRead(U,&u));
33   if (actx->mode == 1) {
34     fvalue[0] = u[1]-actx->lambda1*u[0];
35   }else if (actx->mode == 2) {
36     fvalue[0] = u[1]-actx->lambda2*u[0];
37   }
38   PetscCall(VecRestoreArrayRead(U,&u));
39   PetscFunctionReturn(0);
40 }
41 
42 PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
43 {
44   AppCtx         *actx=(AppCtx*)ctx;
45 
46   PetscFunctionBegin;
47   if (actx->mode == 1) {
48     actx->mode = 2;
49     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t));
50   } else if (actx->mode == 2) {
51     actx->mode = 1;
52     PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t));
53   }
54   PetscFunctionReturn(0);
55 }
56 
57 /*
58      Defines the ODE passed to the ODE solver
59 */
60 static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
61 {
62   AppCtx            *actx=(AppCtx*)ctx;
63   PetscScalar       *f;
64   const PetscScalar *u,*udot;
65 
66   PetscFunctionBegin;
67   /*  The next three lines allow us to access the entries of the vectors directly */
68   PetscCall(VecGetArrayRead(U,&u));
69   PetscCall(VecGetArrayRead(Udot,&udot));
70   PetscCall(VecGetArray(F,&f));
71 
72   if (actx->mode == 1) {
73     f[0] = udot[0]-u[0]+100*u[1];
74     f[1] = udot[1]-10*u[0]-u[1];
75   } else if (actx->mode == 2) {
76     f[0] = udot[0]-u[0]-10*u[1];
77     f[1] = udot[1]+100*u[0]-u[1];
78   }
79 
80   PetscCall(VecRestoreArrayRead(U,&u));
81   PetscCall(VecRestoreArrayRead(Udot,&udot));
82   PetscCall(VecRestoreArray(F,&f));
83   PetscFunctionReturn(0);
84 }
85 
86 /*
87      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
88 */
89 static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
90 {
91   AppCtx            *actx=(AppCtx*)ctx;
92   PetscInt          rowcol[] = {0,1};
93   PetscScalar       J[2][2];
94   const PetscScalar *u,*udot;
95 
96   PetscFunctionBegin;
97   PetscCall(VecGetArrayRead(U,&u));
98   PetscCall(VecGetArrayRead(Udot,&udot));
99 
100   if (actx->mode == 1) {
101     J[0][0] = a-1;                       J[0][1] = 100;
102     J[1][0] = -10;                       J[1][1] = a-1;
103   } else if (actx->mode == 2) {
104     J[0][0] = a-1;                       J[0][1] = -10;
105     J[1][0] = 100;                       J[1][1] = a-1;
106   }
107   PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
108 
109   PetscCall(VecRestoreArrayRead(U,&u));
110   PetscCall(VecRestoreArrayRead(Udot,&udot));
111 
112   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
113   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
114   if (A != B) {
115     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
116     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
117   }
118   PetscFunctionReturn(0);
119 }
120 
121 int main(int argc,char **argv)
122 {
123   TS             ts;            /* ODE integrator */
124   Vec            U;             /* solution will be stored here */
125   Mat            A;             /* Jacobian matrix */
126   PetscMPIInt    size;
127   PetscInt       n = 2;
128   PetscScalar    *u;
129   AppCtx         app;
130   PetscInt       direction[1];
131   PetscBool      terminate[1];
132 
133   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134      Initialize program
135      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136   PetscFunctionBeginUser;
137   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
138   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
139   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
140   app.mode = 1;
141   app.lambda1 = 2.75;
142   app.lambda2 = 0.36;
143   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex1 options","");
144   {
145     PetscCall(PetscOptionsReal("-lambda1","","",app.lambda1,&app.lambda1,NULL));
146     PetscCall(PetscOptionsReal("-lambda2","","",app.lambda2,&app.lambda2,NULL));
147   }
148   PetscOptionsEnd();
149 
150   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151     Create necessary matrix and vectors
152     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
154   PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
155   PetscCall(MatSetType(A,MATDENSE));
156   PetscCall(MatSetFromOptions(A));
157   PetscCall(MatSetUp(A));
158 
159   PetscCall(MatCreateVecs(A,&U,NULL));
160 
161   PetscCall(VecGetArray(U,&u));
162   u[0] = 0;
163   u[1] = 1;
164   PetscCall(VecRestoreArray(U,&u));
165   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166      Create timestepping solver context
167      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
169   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
170   PetscCall(TSSetType(ts,TSCN));
171   PetscCall(TSSetIFunction(ts,NULL,(TSIFunction)IFunction,&app));
172   PetscCall(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&app));
173 
174   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175      Set initial conditions
176    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177   PetscCall(TSSetSolution(ts,U));
178 
179   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180      Set solver options
181    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182   PetscCall(TSSetMaxTime(ts,0.125));
183   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
184   PetscCall(TSSetTimeStep(ts,1./256.));
185   PetscCall(TSSetFromOptions(ts));
186 
187   /* Set directions and terminate flags for the two events */
188   direction[0] = 0;
189   terminate[0] = PETSC_FALSE;
190   PetscCall(TSSetEventHandler(ts,1,direction,terminate,EventFunction,PostEventFunction,(void*)&app));
191 
192   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193      Run timestepping solver
194      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195     PetscCall(TSSolve(ts,U));
196 
197   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
199    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200   PetscCall(MatDestroy(&A));
201   PetscCall(VecDestroy(&U));
202   PetscCall(TSDestroy(&ts));
203 
204   PetscCall(PetscFinalize());
205   return(0);
206 }
207 
208 /*TEST
209 
210    build:
211       requires: !complex
212    test:
213       args: -ts_monitor
214 
215    test:
216       suffix: 2
217       args: -ts_monitor_lg_solution -1
218       requires: x
219 
220 TEST*/
221