xref: /petsc/src/ts/tutorials/hybrid/ex1fwd.c (revision 1d27aa22b2f6148b2c4e3f06a75e0638d6493e09)
1c4762a1bSJed Brown static char help[] = "Trajectory sensitivity of a hybrid system with state-dependent switchings.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown   The dynamics is described by the ODE
5c4762a1bSJed Brown                   u_t = A_i u
6c4762a1bSJed Brown 
7c4762a1bSJed Brown   where A_1 = [ 1  -100
8c4762a1bSJed Brown                 10  1  ],
9c4762a1bSJed Brown         A_2 = [ 1    10
10c4762a1bSJed Brown                -100  1 ].
11c4762a1bSJed Brown   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].
12c4762a1bSJed Brown   Initially u=[0 1]^T and i=1.
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   References:
15*1d27aa22SBarry Smith + * - H. Zhang, S. Abhyankar, E. Constantinescu, M. Mihai, Discrete Adjoint Sensitivity Analysis of Hybrid Dynamical Systems With Switching,
16*1d27aa22SBarry Smith       IEEE Transactions on Circuits and Systems I: Regular Papers, 64(5), May 2017
17606c0280SSatish Balay - * - I. A. Hiskens, M.A. Pai, Trajectory Sensitivity Analysis of Hybrid Systems, IEEE Transactions on Circuits and Systems, Vol 47, No 2, February 2000
18c4762a1bSJed Brown */
19c4762a1bSJed Brown 
20c4762a1bSJed Brown #include <petscts.h>
21c4762a1bSJed Brown 
22c4762a1bSJed Brown typedef struct {
23c4762a1bSJed Brown   PetscScalar lambda1;
24c4762a1bSJed Brown   PetscScalar lambda2;
25c4762a1bSJed Brown   PetscInt    mode; /* mode flag*/
26c4762a1bSJed Brown   PetscReal   print_time;
27c4762a1bSJed Brown } AppCtx;
28c4762a1bSJed Brown 
29d71ae5a4SJacob Faibussowitsch PetscErrorCode MyMonitor(TS ts, PetscInt stepnum, PetscReal time, Vec U, void *ctx)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown   AppCtx      *actx = (AppCtx *)ctx;
32c4762a1bSJed Brown   Mat          sp;
33c4762a1bSJed Brown   PetscScalar *u;
34c4762a1bSJed Brown   PetscInt     nump;
35c4762a1bSJed Brown   FILE        *f;
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   PetscFunctionBegin;
38c4762a1bSJed Brown   if (time >= actx->print_time) {
39c4762a1bSJed Brown     actx->print_time += 1. / 256.;
409566063dSJacob Faibussowitsch     PetscCall(TSForwardGetSensitivities(ts, &nump, &sp));
419566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumn(sp, 2, &u));
42c4762a1bSJed Brown     f = fopen("fwd_sp.out", "a");
4363a3b9bcSJacob Faibussowitsch     PetscCall(PetscFPrintf(PETSC_COMM_WORLD, f, "%20.15lf %20.15lf %20.15lf\n", (double)time, (double)PetscRealPart(u[0]), (double)PetscRealPart(u[1])));
449566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumn(sp, &u));
45c4762a1bSJed Brown     fclose(f);
46c4762a1bSJed Brown   }
473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
50d71ae5a4SJacob Faibussowitsch PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx)
51d71ae5a4SJacob Faibussowitsch {
52c4762a1bSJed Brown   AppCtx            *actx = (AppCtx *)ctx;
53c4762a1bSJed Brown   const PetscScalar *u;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   PetscFunctionBegin;
569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
57c4762a1bSJed Brown   if (actx->mode == 1) {
58c4762a1bSJed Brown     fvalue[0] = u[1] - actx->lambda1 * u[0];
59c4762a1bSJed Brown   } else if (actx->mode == 2) {
60c4762a1bSJed Brown     fvalue[0] = u[1] - actx->lambda2 * u[0];
61c4762a1bSJed Brown   }
629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66d71ae5a4SJacob Faibussowitsch PetscErrorCode ShiftGradients(TS ts, Vec U, AppCtx *actx)
67d71ae5a4SJacob Faibussowitsch {
68c4762a1bSJed Brown   Mat          sp;
69c4762a1bSJed Brown   PetscScalar *x;
70c4762a1bSJed Brown   PetscScalar *u;
71c4762a1bSJed Brown   PetscScalar  tmp[2], A1[2][2], A2[2], denorm;
72c4762a1bSJed Brown   PetscInt     nump;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   PetscFunctionBegin;
759566063dSJacob Faibussowitsch   PetscCall(TSForwardGetSensitivities(ts, &nump, &sp));
769566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   if (actx->mode == 1) {
79c4762a1bSJed Brown     denorm   = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]);
80c4762a1bSJed Brown     A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm + 1.;
81c4762a1bSJed Brown     A1[1][0] = -110. * u[0] * (-actx->lambda1) / denorm;
82c4762a1bSJed Brown     A1[0][1] = 110. * u[1] * 1. / denorm;
83c4762a1bSJed Brown     A1[1][1] = -110. * u[0] * 1. / denorm + 1.;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown     A2[0] = 110. * u[1] * (-u[0]) / denorm;
86c4762a1bSJed Brown     A2[1] = -110. * u[0] * (-u[0]) / denorm;
87c4762a1bSJed Brown   } else {
88c4762a1bSJed Brown     denorm   = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]);
89c4762a1bSJed Brown     A1[0][0] = 110. * u[1] * (actx->lambda2) / denorm + 1;
90c4762a1bSJed Brown     A1[1][0] = -110. * u[0] * (actx->lambda2) / denorm;
91c4762a1bSJed Brown     A1[0][1] = -110. * u[1] * 1. / denorm;
92c4762a1bSJed Brown     A1[1][1] = 110. * u[0] * 1. / denorm + 1.;
93c4762a1bSJed Brown 
94c4762a1bSJed Brown     A2[0] = 0;
95c4762a1bSJed Brown     A2[1] = 0;
96c4762a1bSJed Brown   }
97c4762a1bSJed Brown 
989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
99c4762a1bSJed Brown 
1009566063dSJacob Faibussowitsch   PetscCall(MatDenseGetColumn(sp, 0, &x));
101c4762a1bSJed Brown   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
102c4762a1bSJed Brown   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
103c4762a1bSJed Brown   x[0]   = tmp[0];
104c4762a1bSJed Brown   x[1]   = tmp[1];
1059566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreColumn(sp, &x));
106c4762a1bSJed Brown 
1079566063dSJacob Faibussowitsch   PetscCall(MatDenseGetColumn(sp, 1, &x));
108c4762a1bSJed Brown   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
109c4762a1bSJed Brown   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
110c4762a1bSJed Brown   x[0]   = tmp[0];
111c4762a1bSJed Brown   x[1]   = tmp[1];
1129566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreColumn(sp, &x));
113c4762a1bSJed Brown 
1149566063dSJacob Faibussowitsch   PetscCall(MatDenseGetColumn(sp, 2, &x));
115c4762a1bSJed Brown   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
116c4762a1bSJed Brown   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
117c4762a1bSJed Brown   x[0]   = tmp[0] + A2[0];
118c4762a1bSJed Brown   x[1]   = tmp[1] + A2[1];
1199566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreColumn(sp, &x));
120c4762a1bSJed Brown 
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122c4762a1bSJed Brown }
123c4762a1bSJed Brown 
124d71ae5a4SJacob Faibussowitsch PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx)
125d71ae5a4SJacob Faibussowitsch {
126c4762a1bSJed Brown   AppCtx *actx = (AppCtx *)ctx;
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   PetscFunctionBegin;
1299566063dSJacob Faibussowitsch   /* PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); */
1309566063dSJacob Faibussowitsch   PetscCall(ShiftGradients(ts, U, actx));
131c4762a1bSJed Brown   if (actx->mode == 1) {
132c4762a1bSJed Brown     actx->mode = 2;
1339566063dSJacob Faibussowitsch     /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t)); */
134c4762a1bSJed Brown   } else if (actx->mode == 2) {
135c4762a1bSJed Brown     actx->mode = 1;
1369566063dSJacob Faibussowitsch     /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t)); */
137c4762a1bSJed Brown   }
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown /*
142c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
143c4762a1bSJed Brown */
144d71ae5a4SJacob Faibussowitsch static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx)
145d71ae5a4SJacob Faibussowitsch {
146c4762a1bSJed Brown   AppCtx            *actx = (AppCtx *)ctx;
147c4762a1bSJed Brown   PetscScalar       *f;
148c4762a1bSJed Brown   const PetscScalar *u, *udot;
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   PetscFunctionBegin;
151c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1529566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1549566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   if (actx->mode == 1) {
157c4762a1bSJed Brown     f[0] = udot[0] - u[0] + 100 * u[1];
158c4762a1bSJed Brown     f[1] = udot[1] - 10 * u[0] - u[1];
159c4762a1bSJed Brown   } else if (actx->mode == 2) {
160c4762a1bSJed Brown     f[0] = udot[0] - u[0] - 10 * u[1];
161c4762a1bSJed Brown     f[1] = udot[1] + 100 * u[0] - u[1];
162c4762a1bSJed Brown   }
163c4762a1bSJed Brown 
1649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown /*
171c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
172c4762a1bSJed Brown */
173d71ae5a4SJacob Faibussowitsch static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx)
174d71ae5a4SJacob Faibussowitsch {
175c4762a1bSJed Brown   AppCtx            *actx     = (AppCtx *)ctx;
176c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
177c4762a1bSJed Brown   PetscScalar        J[2][2];
178c4762a1bSJed Brown   const PetscScalar *u, *udot;
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   PetscFunctionBegin;
1819566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1829566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   if (actx->mode == 1) {
1859371c9d4SSatish Balay     J[0][0] = a - 1;
1869371c9d4SSatish Balay     J[0][1] = 100;
1879371c9d4SSatish Balay     J[1][0] = -10;
1889371c9d4SSatish Balay     J[1][1] = a - 1;
189c4762a1bSJed Brown   } else if (actx->mode == 2) {
1909371c9d4SSatish Balay     J[0][0] = a - 1;
1919371c9d4SSatish Balay     J[0][1] = -10;
1929371c9d4SSatish Balay     J[1][0] = 100;
1939371c9d4SSatish Balay     J[1][1] = a - 1;
194c4762a1bSJed Brown   }
1959566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
196c4762a1bSJed Brown 
1979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
199c4762a1bSJed Brown 
2009566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2019566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
202c4762a1bSJed Brown   if (A != B) {
2039566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2049566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
205c4762a1bSJed Brown   }
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
207c4762a1bSJed Brown }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown /* Matrix JacobianP is constant so that it only needs to be evaluated once */
210d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat Ap, void *ctx)
211d71ae5a4SJacob Faibussowitsch {
212c4762a1bSJed Brown   PetscFunctionBeginUser;
2133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
214c4762a1bSJed Brown }
215c4762a1bSJed Brown 
216d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
217d71ae5a4SJacob Faibussowitsch {
218c4762a1bSJed Brown   TS           ts; /* ODE integrator */
219c4762a1bSJed Brown   Vec          U;  /* solution will be stored here */
220c4762a1bSJed Brown   Mat          A;  /* Jacobian matrix */
221c4762a1bSJed Brown   Mat          Ap; /* Jacobian dfdp */
222c4762a1bSJed Brown   PetscMPIInt  size;
223c4762a1bSJed Brown   PetscInt     n = 2;
224c4762a1bSJed Brown   PetscScalar *u;
225c4762a1bSJed Brown   AppCtx       app;
226c4762a1bSJed Brown   PetscInt     direction[1];
227c4762a1bSJed Brown   PetscBool    terminate[1];
228c4762a1bSJed Brown   Mat          sp;
229c4762a1bSJed Brown   PetscReal    tend;
230c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
231c4762a1bSJed Brown      Initialize program
232c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
233327415f7SBarry Smith   PetscFunctionBeginUser;
2349566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
2359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2363c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
237c4762a1bSJed Brown   app.mode       = 1;
238c4762a1bSJed Brown   app.lambda1    = 2.75;
239c4762a1bSJed Brown   app.lambda2    = 0.36;
240c4762a1bSJed Brown   app.print_time = 1. / 256.;
241c4762a1bSJed Brown   tend           = 0.125;
242d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1fwd options", "");
243c4762a1bSJed Brown   {
2449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL));
2459566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL));
2469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tend", "", "", tend, &tend, NULL));
247c4762a1bSJed Brown   }
248d0609cedSBarry Smith   PetscOptionsEnd();
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
251c4762a1bSJed Brown     Create necessary matrix and vectors
252c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2539566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2549566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
2559566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATDENSE));
2569566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2579566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
258c4762a1bSJed Brown 
2599566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &U, NULL));
260c4762a1bSJed Brown 
2619566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap));
2629566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Ap, n, 3, PETSC_DETERMINE, PETSC_DETERMINE));
2639566063dSJacob Faibussowitsch   PetscCall(MatSetType(Ap, MATDENSE));
2649566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Ap));
2659566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Ap));
2669566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(Ap));
267c4762a1bSJed Brown 
2689566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, n, 3, NULL, &sp));
2699566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(sp));
2709566063dSJacob Faibussowitsch   PetscCall(MatShift(sp, 1.0));
271c4762a1bSJed Brown 
2729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
273c4762a1bSJed Brown   u[0] = 0;
274c4762a1bSJed Brown   u[1] = 1;
2759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
278c4762a1bSJed Brown      Create timestepping solver context
279c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2809566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2819566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
2829566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSCN));
2839566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &app));
2849566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &app));
285c4762a1bSJed Brown 
286c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
287c4762a1bSJed Brown      Set initial conditions
288c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2899566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
290c4762a1bSJed Brown 
2919566063dSJacob Faibussowitsch   PetscCall(TSForwardSetSensitivities(ts, 3, sp));
292c4762a1bSJed Brown   /*   Set RHS JacobianP */
2939566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ts, Ap, RHSJacobianP, &app));
294c4762a1bSJed Brown 
295c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
296c4762a1bSJed Brown      Set solver options
297c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2989566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, tend));
2999566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
3009566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 1. / 256.));
301f3fa974cSJacob Faibussowitsch   PetscCall(TSMonitorSet(ts, MyMonitor, &app, NULL));
3029566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   /* Set directions and terminate flags for the two events */
305c4762a1bSJed Brown   direction[0] = 0;
306c4762a1bSJed Brown   terminate[0] = PETSC_FALSE;
3079566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app));
308c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
309c4762a1bSJed Brown      Run timestepping solver
310c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3119566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
314c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
315c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3189566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
319c4762a1bSJed Brown 
3209566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ap));
3219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&sp));
3229566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
323d0609cedSBarry Smith   return (0);
324c4762a1bSJed Brown }
325c4762a1bSJed Brown 
326c4762a1bSJed Brown /*TEST
327c4762a1bSJed Brown 
328c4762a1bSJed Brown    build:
329c4762a1bSJed Brown       requires: !complex
330c4762a1bSJed Brown 
331c4762a1bSJed Brown    test:
332c4762a1bSJed Brown       args: -ts_monitor
333c4762a1bSJed Brown 
334c4762a1bSJed Brown TEST*/
335