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