xref: /petsc/src/ts/tutorials/hybrid/ex1adj.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown static char help[] = "Adjoint 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:
15606c0280SSatish Balay + * - 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
16606c0280SSatish 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
17c4762a1bSJed Brown */
18c4762a1bSJed Brown 
19c4762a1bSJed Brown #include <petscts.h>
20c4762a1bSJed Brown 
21c4762a1bSJed Brown typedef struct {
22c4762a1bSJed Brown   PetscScalar lambda1;
23c4762a1bSJed Brown   PetscScalar lambda2;
24c4762a1bSJed Brown   PetscInt    mode; /* mode flag*/
25c4762a1bSJed Brown } AppCtx;
26c4762a1bSJed Brown 
279371c9d4SSatish Balay PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscScalar *fvalue, void *ctx) {
28c4762a1bSJed Brown   AppCtx            *actx = (AppCtx *)ctx;
29c4762a1bSJed Brown   const PetscScalar *u;
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
33c4762a1bSJed Brown   if (actx->mode == 1) {
34c4762a1bSJed Brown     fvalue[0] = u[1] - actx->lambda1 * u[0];
35c4762a1bSJed Brown   } else if (actx->mode == 2) {
36c4762a1bSJed Brown     fvalue[0] = u[1] - actx->lambda2 * u[0];
37c4762a1bSJed Brown   }
389566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
39c4762a1bSJed Brown   PetscFunctionReturn(0);
40c4762a1bSJed Brown }
41c4762a1bSJed Brown 
429371c9d4SSatish Balay PetscErrorCode ShiftGradients(TS ts, Vec U, AppCtx *actx) {
43c4762a1bSJed Brown   Vec               *lambda, *mu;
44c4762a1bSJed Brown   PetscScalar       *x, *y;
45c4762a1bSJed Brown   const PetscScalar *u;
46c4762a1bSJed Brown   PetscScalar        tmp[2], A1[2][2], A2[2], denorm;
47c4762a1bSJed Brown   PetscInt           numcost;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   PetscFunctionBegin;
509566063dSJacob Faibussowitsch   PetscCall(TSGetCostGradients(ts, &numcost, &lambda, &mu));
519566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   if (actx->mode == 2) {
54c4762a1bSJed Brown     denorm   = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]);
55c4762a1bSJed Brown     A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm + 1.;
56c4762a1bSJed Brown     A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm;
57c4762a1bSJed Brown     A1[1][0] = 110. * u[1] * 1. / denorm;
58c4762a1bSJed Brown     A1[1][1] = -110. * u[0] * 1. / denorm + 1.;
59c4762a1bSJed Brown 
60c4762a1bSJed Brown     A2[0] = 110. * u[1] * (-u[0]) / denorm;
61c4762a1bSJed Brown     A2[1] = -110. * u[0] * (-u[0]) / denorm;
62c4762a1bSJed Brown   } else {
63c4762a1bSJed Brown     denorm   = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]);
64c4762a1bSJed Brown     A1[0][0] = 110. * u[1] * (actx->lambda2) / denorm + 1;
65c4762a1bSJed Brown     A1[0][1] = -110. * u[0] * (actx->lambda2) / denorm;
66c4762a1bSJed Brown     A1[1][0] = -110. * u[1] * 1. / denorm;
67c4762a1bSJed Brown     A1[1][1] = 110. * u[0] * 1. / denorm + 1.;
68c4762a1bSJed Brown 
69c4762a1bSJed Brown     A2[0] = 0;
70c4762a1bSJed Brown     A2[1] = 0;
71c4762a1bSJed Brown   }
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0], &x));
769566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0], &y));
77c4762a1bSJed Brown   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
78c4762a1bSJed Brown   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
79c4762a1bSJed Brown   y[0]   = y[0] + A2[0] * x[0] + A2[1] * x[1];
80c4762a1bSJed Brown   x[0]   = tmp[0];
81c4762a1bSJed Brown   x[1]   = tmp[1];
829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0], &y));
839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0], &x));
84c4762a1bSJed Brown 
859566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[1], &x));
869566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[1], &y));
87c4762a1bSJed Brown   tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
88c4762a1bSJed Brown   tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
89c4762a1bSJed Brown   y[0]   = y[0] + A2[0] * x[0] + A2[1] * x[1];
90c4762a1bSJed Brown   x[0]   = tmp[0];
91c4762a1bSJed Brown   x[1]   = tmp[1];
929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[1], &y));
939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[1], &x));
94c4762a1bSJed Brown   PetscFunctionReturn(0);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
979371c9d4SSatish Balay PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx) {
98c4762a1bSJed Brown   AppCtx *actx = (AppCtx *)ctx;
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   PetscFunctionBegin;
1019566063dSJacob Faibussowitsch   /* PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); */
102*48a46eb9SPierre Jolivet   if (!forwardsolve) PetscCall(ShiftGradients(ts, U, actx));
103c4762a1bSJed Brown   if (actx->mode == 1) {
104c4762a1bSJed Brown     actx->mode = 2;
1059566063dSJacob Faibussowitsch     /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t)); */
106c4762a1bSJed Brown   } else if (actx->mode == 2) {
107c4762a1bSJed Brown     actx->mode = 1;
1089566063dSJacob Faibussowitsch     /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t)); */
109c4762a1bSJed Brown   }
110c4762a1bSJed Brown   PetscFunctionReturn(0);
111c4762a1bSJed Brown }
112c4762a1bSJed Brown 
113c4762a1bSJed Brown /*
114c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
115c4762a1bSJed Brown */
1169371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx) {
117c4762a1bSJed Brown   AppCtx            *actx = (AppCtx *)ctx;
118c4762a1bSJed Brown   PetscScalar       *f;
119c4762a1bSJed Brown   const PetscScalar *u, *udot;
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   PetscFunctionBegin;
122c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1239566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1259566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   if (actx->mode == 1) {
128c4762a1bSJed Brown     f[0] = udot[0] - u[0] + 100 * u[1];
129c4762a1bSJed Brown     f[1] = udot[1] - 10 * u[0] - u[1];
130c4762a1bSJed Brown   } else if (actx->mode == 2) {
131c4762a1bSJed Brown     f[0] = udot[0] - u[0] - 10 * u[1];
132c4762a1bSJed Brown     f[1] = udot[1] + 100 * u[0] - u[1];
133c4762a1bSJed Brown   }
134c4762a1bSJed Brown 
1359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
138c4762a1bSJed Brown   PetscFunctionReturn(0);
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown /*
142c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
143c4762a1bSJed Brown */
1449371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) {
145c4762a1bSJed Brown   AppCtx            *actx     = (AppCtx *)ctx;
146c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
147c4762a1bSJed Brown   PetscScalar        J[2][2];
148c4762a1bSJed Brown   const PetscScalar *u, *udot;
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   PetscFunctionBegin;
1519566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1529566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   if (actx->mode == 1) {
1559371c9d4SSatish Balay     J[0][0] = a - 1;
1569371c9d4SSatish Balay     J[0][1] = 100;
1579371c9d4SSatish Balay     J[1][0] = -10;
1589371c9d4SSatish Balay     J[1][1] = a - 1;
159c4762a1bSJed Brown   } else if (actx->mode == 2) {
1609371c9d4SSatish Balay     J[0][0] = a - 1;
1619371c9d4SSatish Balay     J[0][1] = -10;
1629371c9d4SSatish Balay     J[1][0] = 100;
1639371c9d4SSatish Balay     J[1][1] = a - 1;
164c4762a1bSJed Brown   }
1659566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
166c4762a1bSJed Brown 
1679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
169c4762a1bSJed Brown 
1709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
172c4762a1bSJed Brown   if (A != B) {
1739566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1749566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
175c4762a1bSJed Brown   }
176c4762a1bSJed Brown   PetscFunctionReturn(0);
177c4762a1bSJed Brown }
178c4762a1bSJed Brown 
179c4762a1bSJed Brown /* Matrix JacobianP is constant so that it only needs to be evaluated once */
1809371c9d4SSatish Balay static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx) {
181c4762a1bSJed Brown   PetscFunctionBeginUser;
182c4762a1bSJed Brown   PetscFunctionReturn(0);
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
1859371c9d4SSatish Balay int main(int argc, char **argv) {
186c4762a1bSJed Brown   TS           ts; /* ODE integrator */
187c4762a1bSJed Brown   Vec          U;  /* solution will be stored here */
188c4762a1bSJed Brown   Mat          A;  /* Jacobian matrix */
189c4762a1bSJed Brown   Mat          Ap; /* dfdp */
190c4762a1bSJed Brown   PetscMPIInt  size;
191c4762a1bSJed Brown   PetscInt     n = 2;
192c4762a1bSJed Brown   PetscScalar *u, *v;
193c4762a1bSJed Brown   AppCtx       app;
194c4762a1bSJed Brown   PetscInt     direction[1];
195c4762a1bSJed Brown   PetscBool    terminate[1];
196c4762a1bSJed Brown   Vec          lambda[2], mu[2];
197c4762a1bSJed Brown   PetscReal    tend;
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   FILE *f;
200c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201c4762a1bSJed Brown      Initialize program
202c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203327415f7SBarry Smith   PetscFunctionBeginUser;
2049566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
2059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2063c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
207c4762a1bSJed Brown   app.mode    = 1;
208c4762a1bSJed Brown   app.lambda1 = 2.75;
209c4762a1bSJed Brown   app.lambda2 = 0.36;
210c4762a1bSJed Brown   tend        = 0.125;
211d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1adj options", "");
212c4762a1bSJed Brown   {
2139566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL));
2149566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL));
2159566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tend", "", "", tend, &tend, NULL));
216c4762a1bSJed Brown   }
217d0609cedSBarry Smith   PetscOptionsEnd();
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220c4762a1bSJed Brown     Create necessary matrix and vectors
221c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2229566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2239566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
2249566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATDENSE));
2259566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2269566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
227c4762a1bSJed Brown 
2289566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &U, NULL));
229c4762a1bSJed Brown 
2309566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap));
2319566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Ap, n, 1, PETSC_DETERMINE, PETSC_DETERMINE));
2329566063dSJacob Faibussowitsch   PetscCall(MatSetType(Ap, MATDENSE));
2339566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Ap));
2349566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Ap));
2359566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(Ap)); /* initialize to zeros */
236c4762a1bSJed Brown 
2379566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
238c4762a1bSJed Brown   u[0] = 0;
239c4762a1bSJed Brown   u[1] = 1;
2409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
241c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
242c4762a1bSJed Brown      Create timestepping solver context
243c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2449566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2459566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
2469566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSCN));
2479566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &app));
2489566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &app));
2499566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ts, Ap, RHSJacobianP, &app));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252c4762a1bSJed Brown      Set initial conditions
253c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2549566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
257c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
258c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2599566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262c4762a1bSJed Brown      Set solver options
263c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2649566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, tend));
2659566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
2669566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 1. / 256.));
2679566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /* Set directions and terminate flags for the two events */
270c4762a1bSJed Brown   direction[0] = 0;
271c4762a1bSJed Brown   terminate[0] = PETSC_FALSE;
2729566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app));
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
275c4762a1bSJed Brown      Run timestepping solver
276c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2779566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
280c4762a1bSJed Brown      Adjoint model starts here
281c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2829566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &lambda[0], NULL));
2839566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &lambda[1], NULL));
284c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
2859566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(lambda[0]));
2869566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(lambda[1]));
2879566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0], &u));
288c4762a1bSJed Brown   u[0] = 1.;
2899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0], &u));
2909566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[1], &u));
291c4762a1bSJed Brown   u[1] = 1.;
2929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[1], &u));
293c4762a1bSJed Brown 
2949566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Ap, &mu[0], NULL));
2959566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Ap, &mu[1], NULL));
2969566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(mu[0]));
2979566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(mu[1]));
2989566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 2, lambda, mu));
299c4762a1bSJed Brown 
3009566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   /*
3039566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD));
3049566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD));
3059566063dSJacob Faibussowitsch   PetscCall(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD));
3069566063dSJacob Faibussowitsch   PetscCall(VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD));
307c4762a1bSJed Brown   */
3089566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0], &u));
3099566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[1], &v));
310c4762a1bSJed Brown   f = fopen("adj_mu.out", "a");
31163a3b9bcSJacob Faibussowitsch   PetscCall(PetscFPrintf(PETSC_COMM_WORLD, f, "%20.15lf %20.15lf %20.15lf\n", (double)tend, (double)PetscRealPart(u[0]), (double)PetscRealPart(v[0])));
3129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0], &u));
3139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[1], &v));
314c4762a1bSJed Brown   fclose(f);
315c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
316c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
317c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3209566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
321c4762a1bSJed Brown 
3229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ap));
3239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[0]));
3249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[1]));
3259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[0]));
3269566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[1]));
3279566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
328b122ec5aSJacob Faibussowitsch   return 0;
329c4762a1bSJed Brown }
330c4762a1bSJed Brown 
331c4762a1bSJed Brown /*TEST
332c4762a1bSJed Brown 
333c4762a1bSJed Brown    build:
334c4762a1bSJed Brown       requires: !complex
335c4762a1bSJed Brown 
336c4762a1bSJed Brown    test:
337c4762a1bSJed Brown       args: -ts_monitor -ts_adjoint_monitor
338c4762a1bSJed Brown 
339c4762a1bSJed Brown TEST*/
340