xref: /petsc/src/ts/tutorials/hybrid/ex1fd.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown static char help[] = "An example of hybrid system using TS event.\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   Reference:
15c4762a1bSJed Brown   I. A. Hiskens, M.A. Pai, Trajectory Sensitivity Analysis of Hybrid Systems, IEEE Transactions on Circuits and Systems, Vol 47, No 2, February 2000
16c4762a1bSJed Brown */
17c4762a1bSJed Brown 
18c4762a1bSJed Brown #include <petscts.h>
19c4762a1bSJed Brown 
20c4762a1bSJed Brown typedef struct {
21c4762a1bSJed Brown   PetscReal lambda1;
22c4762a1bSJed Brown   PetscReal lambda2;
23c4762a1bSJed Brown   PetscInt  mode; /* mode flag*/
24c4762a1bSJed Brown } AppCtx;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown PetscErrorCode FWDRun(TS, Vec, void *);
27c4762a1bSJed Brown 
28*d71ae5a4SJacob Faibussowitsch PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx)
29*d71ae5a4SJacob Faibussowitsch {
30c4762a1bSJed Brown   AppCtx            *actx = (AppCtx *)ctx;
31c4762a1bSJed Brown   const PetscScalar *u;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   PetscFunctionBegin;
349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
35c4762a1bSJed Brown   if (actx->mode == 1) {
36c4762a1bSJed Brown     fvalue[0] = u[1] - actx->lambda1 * u[0];
37c4762a1bSJed Brown   } else if (actx->mode == 2) {
38c4762a1bSJed Brown     fvalue[0] = u[1] - actx->lambda2 * u[0];
39c4762a1bSJed Brown   }
409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
41c4762a1bSJed Brown   PetscFunctionReturn(0);
42c4762a1bSJed Brown }
43c4762a1bSJed Brown 
44*d71ae5a4SJacob Faibussowitsch PetscErrorCode ShiftGradients(TS ts, Vec U, AppCtx *actx)
45*d71ae5a4SJacob Faibussowitsch {
46c4762a1bSJed Brown   Vec               *lambda, *mu;
47c4762a1bSJed Brown   PetscScalar       *x, *y;
48c4762a1bSJed Brown   const PetscScalar *u;
49c4762a1bSJed Brown   PetscScalar        tmp[2], A1[2][2], A2[2], denorm1, denorm2;
50c4762a1bSJed Brown   PetscInt           numcost;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   PetscFunctionBegin;
539566063dSJacob Faibussowitsch   PetscCall(TSGetCostGradients(ts, &numcost, &lambda, &mu));
549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   if (actx->mode == 2) {
57c4762a1bSJed Brown     denorm1  = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]);
58c4762a1bSJed Brown     denorm2  = -actx->lambda1 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]);
59c4762a1bSJed Brown     A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm1 + 1.;
60c4762a1bSJed Brown     A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm1;
61c4762a1bSJed Brown     A1[1][0] = 110. * u[1] * 1. / denorm1;
62c4762a1bSJed Brown     A1[1][1] = -110. * u[0] * 1. / denorm1 + 1.;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown     A2[0] = 110. * u[1] * (-u[0]) / denorm2;
65c4762a1bSJed Brown     A2[1] = -110. * u[0] * (-u[0]) / denorm2;
66c4762a1bSJed Brown   } else {
67c4762a1bSJed Brown     denorm2  = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]);
68c4762a1bSJed Brown     A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm2 + 1;
69c4762a1bSJed Brown     A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm2;
70c4762a1bSJed Brown     A1[1][0] = 110. * u[1] * 1. / denorm2;
71c4762a1bSJed Brown     A1[1][1] = -110. * u[0] * 1. / denorm2 + 1.;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown     A2[0] = 0;
74c4762a1bSJed Brown     A2[1] = 0;
75c4762a1bSJed Brown   }
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0], &x));
809566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0], &y));
81c4762a1bSJed Brown   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
82c4762a1bSJed Brown   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
83c4762a1bSJed Brown   y[0]   = y[0] + A2[0] * x[0] + A2[1] * x[1];
84c4762a1bSJed Brown   x[0]   = tmp[0];
85c4762a1bSJed Brown   x[1]   = tmp[1];
869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0], &y));
879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0], &x));
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[1], &x));
909566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[1], &y));
91c4762a1bSJed Brown   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
92c4762a1bSJed Brown   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
93c4762a1bSJed Brown   y[0]   = y[0] + A2[0] * x[0] + A2[1] * x[1];
94c4762a1bSJed Brown   x[0]   = tmp[0];
95c4762a1bSJed Brown   x[1]   = tmp[1];
969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[1], &y));
979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[1], &x));
98c4762a1bSJed Brown   PetscFunctionReturn(0);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101*d71ae5a4SJacob Faibussowitsch PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx)
102*d71ae5a4SJacob Faibussowitsch {
103c4762a1bSJed Brown   AppCtx *actx = (AppCtx *)ctx;
104c4762a1bSJed Brown 
105c4762a1bSJed Brown   PetscFunctionBegin;
10648a46eb9SPierre Jolivet   if (!forwardsolve) PetscCall(ShiftGradients(ts, U, actx));
107c4762a1bSJed Brown   if (actx->mode == 1) {
108c4762a1bSJed Brown     actx->mode = 2;
1099566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Change from mode 1 to 2 at t = %f \n", (double)t));
110c4762a1bSJed Brown   } else if (actx->mode == 2) {
111c4762a1bSJed Brown     actx->mode = 1;
1129566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Change from mode 2 to 1 at t = %f \n", (double)t));
113c4762a1bSJed Brown   }
114c4762a1bSJed Brown   PetscFunctionReturn(0);
115c4762a1bSJed Brown }
116c4762a1bSJed Brown 
117c4762a1bSJed Brown /*
118c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
119c4762a1bSJed Brown */
120*d71ae5a4SJacob Faibussowitsch static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx)
121*d71ae5a4SJacob Faibussowitsch {
122c4762a1bSJed Brown   AppCtx            *actx = (AppCtx *)ctx;
123c4762a1bSJed Brown   PetscScalar       *f;
124c4762a1bSJed Brown   const PetscScalar *u, *udot;
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   PetscFunctionBegin;
127c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1289566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1309566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   if (actx->mode == 1) {
133c4762a1bSJed Brown     f[0] = udot[0] - u[0] + 100 * u[1];
134c4762a1bSJed Brown     f[1] = udot[1] - 10 * u[0] - u[1];
135c4762a1bSJed Brown   } else if (actx->mode == 2) {
136c4762a1bSJed Brown     f[0] = udot[0] - u[0] - 10 * u[1];
137c4762a1bSJed Brown     f[1] = udot[1] + 100 * u[0] - u[1];
138c4762a1bSJed Brown   }
139c4762a1bSJed Brown 
1409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
143c4762a1bSJed Brown   PetscFunctionReturn(0);
144c4762a1bSJed Brown }
145c4762a1bSJed Brown 
146c4762a1bSJed Brown /*
147c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
148c4762a1bSJed Brown */
149*d71ae5a4SJacob Faibussowitsch static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx)
150*d71ae5a4SJacob Faibussowitsch {
151c4762a1bSJed Brown   AppCtx            *actx     = (AppCtx *)ctx;
152c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
153c4762a1bSJed Brown   PetscScalar        J[2][2];
154c4762a1bSJed Brown   const PetscScalar *u, *udot;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   PetscFunctionBegin;
1579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   if (actx->mode == 1) {
1619371c9d4SSatish Balay     J[0][0] = a - 1;
1629371c9d4SSatish Balay     J[0][1] = 100;
1639371c9d4SSatish Balay     J[1][0] = -10;
1649371c9d4SSatish Balay     J[1][1] = a - 1;
165c4762a1bSJed Brown   } else if (actx->mode == 2) {
1669371c9d4SSatish Balay     J[0][0] = a - 1;
1679371c9d4SSatish Balay     J[0][1] = -10;
1689371c9d4SSatish Balay     J[1][0] = 100;
1699371c9d4SSatish Balay     J[1][1] = a - 1;
170c4762a1bSJed Brown   }
1719566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
172c4762a1bSJed Brown 
1739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
175c4762a1bSJed Brown 
1769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1779566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
178c4762a1bSJed Brown   if (A != B) {
1799566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1809566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
181c4762a1bSJed Brown   }
182c4762a1bSJed Brown   PetscFunctionReturn(0);
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
185*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
186*d71ae5a4SJacob Faibussowitsch {
187c4762a1bSJed Brown   TS           ts; /* ODE integrator */
188c4762a1bSJed Brown   Vec          U;  /* solution will be stored here */
189c4762a1bSJed Brown   Mat          A;  /* Jacobian matrix */
190c4762a1bSJed Brown   Mat          Ap; /* dfdp */
191c4762a1bSJed Brown   PetscMPIInt  size;
192c4762a1bSJed Brown   PetscInt     n = 2;
193c4762a1bSJed Brown   PetscScalar *u;
194c4762a1bSJed Brown   AppCtx       app;
195c4762a1bSJed Brown   PetscInt     direction[1];
196c4762a1bSJed Brown   PetscBool    terminate[1];
197c4762a1bSJed Brown   PetscReal    delta, tmp[2], sensi[2];
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   delta = 1e-8;
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202c4762a1bSJed Brown      Initialize program
203c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
204327415f7SBarry Smith   PetscFunctionBeginUser;
2059566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
2069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2073c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
208c4762a1bSJed Brown   app.mode    = 1;
209c4762a1bSJed Brown   app.lambda1 = 2.75;
210c4762a1bSJed Brown   app.lambda2 = 0.36;
211d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1 options", "");
212c4762a1bSJed Brown   {
2139566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL));
2149566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL));
215c4762a1bSJed Brown   }
216d0609cedSBarry Smith   PetscOptionsEnd();
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219c4762a1bSJed Brown     Create necessary matrix and vectors
220c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2219566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2229566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
2239566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATDENSE));
2249566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2259566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
226c4762a1bSJed Brown 
2279566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &U, NULL));
228c4762a1bSJed Brown 
2299566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap));
2309566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Ap, n, 1, PETSC_DETERMINE, PETSC_DETERMINE));
2319566063dSJacob Faibussowitsch   PetscCall(MatSetType(Ap, MATDENSE));
2329566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Ap));
2339566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Ap));
2349566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(Ap)); /* initialize to zeros */
235c4762a1bSJed Brown 
2369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
237c4762a1bSJed Brown   u[0] = 0;
238c4762a1bSJed Brown   u[1] = 1;
2399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
240c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
241c4762a1bSJed Brown      Create timestepping solver context
242c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2439566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2449566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
2459566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSCN));
2469566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &app));
2479566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &app));
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
250c4762a1bSJed Brown      Set initial conditions
251c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2529566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
255c4762a1bSJed Brown      Set solver options
256c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2579566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 0.125));
2589566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
2599566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 1. / 256.));
2609566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /* Set directions and terminate flags for the two events */
263c4762a1bSJed Brown   direction[0] = 0;
264c4762a1bSJed Brown   terminate[0] = PETSC_FALSE;
2659566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app));
266c4762a1bSJed Brown 
267c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
268c4762a1bSJed Brown      Run timestepping solver
269c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2709566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
271c4762a1bSJed Brown 
2729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
273c4762a1bSJed Brown   tmp[0] = u[0];
274c4762a1bSJed Brown   tmp[1] = u[1];
275c4762a1bSJed Brown 
276c4762a1bSJed Brown   u[0] = 0 + delta;
277c4762a1bSJed Brown   u[1] = 1;
2789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
279c4762a1bSJed Brown 
2809566063dSJacob Faibussowitsch   PetscCall(FWDRun(ts, U, (void *)&app));
281c4762a1bSJed Brown 
2829566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
283c4762a1bSJed Brown   sensi[0] = (u[0] - tmp[0]) / delta;
284c4762a1bSJed Brown   sensi[1] = (u[1] - tmp[1]) / delta;
28563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "d x1(tf) /d x1(t0) = %f d x2(tf) / d x1(t0) = %f \n", (double)sensi[0], (double)sensi[1]));
286c4762a1bSJed Brown   u[0] = 0;
287c4762a1bSJed Brown   u[1] = 1 + delta;
2889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
289c4762a1bSJed Brown 
2909566063dSJacob Faibussowitsch   PetscCall(FWDRun(ts, U, (void *)&app));
291c4762a1bSJed Brown 
2929566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
293c4762a1bSJed Brown   sensi[0] = (u[0] - tmp[0]) / delta;
294c4762a1bSJed Brown   sensi[1] = (u[1] - tmp[1]) / delta;
29563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "d x1(tf) /d x2(t0) = %f d x2(tf) / d x2(t0) = %f \n", (double)sensi[0], (double)sensi[1]));
296c4762a1bSJed Brown   u[0]        = 0;
297c4762a1bSJed Brown   u[1]        = 1;
298c4762a1bSJed Brown   app.lambda1 = app.lambda1 + delta;
2999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
300c4762a1bSJed Brown 
3019566063dSJacob Faibussowitsch   PetscCall(FWDRun(ts, U, (void *)&app));
302c4762a1bSJed Brown 
3039566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
304c4762a1bSJed Brown   sensi[0] = (u[0] - tmp[0]) / delta;
305c4762a1bSJed Brown   sensi[1] = (u[1] - tmp[1]) / delta;
30663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "Final gradients: d x1(tf) /d p = %f d x2(tf) / d p = %f \n", (double)sensi[0], (double)sensi[1]));
3079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
310c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
311c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3149566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
315c4762a1bSJed Brown 
3169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ap));
3179566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
318b122ec5aSJacob Faibussowitsch   return 0;
319c4762a1bSJed Brown }
320c4762a1bSJed Brown 
321*d71ae5a4SJacob Faibussowitsch PetscErrorCode FWDRun(TS ts, Vec U0, void *ctx0)
322*d71ae5a4SJacob Faibussowitsch {
323c4762a1bSJed Brown   Vec     U; /* solution will be stored here */
324c4762a1bSJed Brown   AppCtx *ctx = (AppCtx *)ctx0;
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   PetscFunctionBeginUser;
3279566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &U));
3289566063dSJacob Faibussowitsch   PetscCall(VecCopy(U0, U));
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   ctx->mode = 1;
331c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
332c4762a1bSJed Brown      Run timestepping solver
333c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3349566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
335c4762a1bSJed Brown 
3369566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
337c4762a1bSJed Brown 
338c4762a1bSJed Brown   PetscFunctionReturn(0);
339c4762a1bSJed Brown }
340c4762a1bSJed Brown 
341c4762a1bSJed Brown /*TEST
342c4762a1bSJed Brown 
343c4762a1bSJed Brown    build:
344dfd57a17SPierre Jolivet       requires: !defined(PETSC_USE_CXXCOMPLEX)
345c4762a1bSJed Brown 
346c4762a1bSJed Brown    test:
347c4762a1bSJed Brown       args: -ts_event_tol 1e-9
348c4762a1bSJed Brown       timeoutfactor: 18
349c4762a1bSJed Brown       requires: !single
350c4762a1bSJed Brown 
351c4762a1bSJed Brown TEST*/
352