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:
151d27aa22SBarry Smith + * - H. Zhang, S. Abhyankar, E. Constantinescu, M. Mihai, Discrete Adjoint Sensitivity Analysis of Hybrid Dynamical Systems With Switching,
161d27aa22SBarry 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 } AppCtx;
27c4762a1bSJed Brown
EventFunction(TS ts,PetscReal t,Vec U,PetscReal * fvalue,PetscCtx ctx)28*2a8381b2SBarry Smith PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal *fvalue, PetscCtx ctx)
29d71ae5a4SJacob 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) {
36fa9584fbSIlya Fursov fvalue[0] = PetscRealPart(u[1] - actx->lambda1 * u[0]);
37c4762a1bSJed Brown } else if (actx->mode == 2) {
38fa9584fbSIlya Fursov fvalue[0] = PetscRealPart(u[1] - actx->lambda2 * u[0]);
39c4762a1bSJed Brown }
409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
42c4762a1bSJed Brown }
43c4762a1bSJed Brown
ShiftGradients(TS ts,Vec U,AppCtx * actx)44d71ae5a4SJacob Faibussowitsch PetscErrorCode ShiftGradients(TS ts, Vec U, AppCtx *actx)
45d71ae5a4SJacob 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], denorm;
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 denorm = -actx->lambda1 * (u[0] - 100. * u[1]) + 1. * (10. * u[0] + u[1]);
58c4762a1bSJed Brown A1[0][0] = 110. * u[1] * (-actx->lambda1) / denorm + 1.;
59c4762a1bSJed Brown A1[0][1] = -110. * u[0] * (-actx->lambda1) / denorm;
60c4762a1bSJed Brown A1[1][0] = 110. * u[1] * 1. / denorm;
61c4762a1bSJed Brown A1[1][1] = -110. * u[0] * 1. / denorm + 1.;
62c4762a1bSJed Brown
63c4762a1bSJed Brown A2[0] = 110. * u[1] * (-u[0]) / denorm;
64c4762a1bSJed Brown A2[1] = -110. * u[0] * (-u[0]) / denorm;
65c4762a1bSJed Brown } else {
66c4762a1bSJed Brown denorm = -actx->lambda2 * (u[0] + 10. * u[1]) + 1. * (-100. * u[0] + u[1]);
67c4762a1bSJed Brown A1[0][0] = 110. * u[1] * (actx->lambda2) / denorm + 1;
68c4762a1bSJed Brown A1[0][1] = -110. * u[0] * (actx->lambda2) / denorm;
69c4762a1bSJed Brown A1[1][0] = -110. * u[1] * 1. / denorm;
70c4762a1bSJed Brown A1[1][1] = 110. * u[0] * 1. / denorm + 1.;
71c4762a1bSJed Brown
72c4762a1bSJed Brown A2[0] = 0;
73c4762a1bSJed Brown A2[1] = 0;
74c4762a1bSJed Brown }
75c4762a1bSJed Brown
769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
77c4762a1bSJed Brown
789566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &x));
799566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &y));
80c4762a1bSJed Brown tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
81c4762a1bSJed Brown tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
82c4762a1bSJed Brown y[0] = y[0] + A2[0] * x[0] + A2[1] * x[1];
83c4762a1bSJed Brown x[0] = tmp[0];
84c4762a1bSJed Brown x[1] = tmp[1];
859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &y));
869566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &x));
87c4762a1bSJed Brown
889566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[1], &x));
899566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[1], &y));
90c4762a1bSJed Brown tmp[0] = A1[0][0] * x[0] + A1[0][1] * x[1];
91c4762a1bSJed Brown tmp[1] = A1[1][0] * x[0] + A1[1][1] * x[1];
92c4762a1bSJed Brown y[0] = y[0] + A2[0] * x[0] + A2[1] * x[1];
93c4762a1bSJed Brown x[0] = tmp[0];
94c4762a1bSJed Brown x[1] = tmp[1];
959566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[1], &y));
969566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[1], &x));
973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
98c4762a1bSJed Brown }
99c4762a1bSJed Brown
PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,PetscCtx ctx)100*2a8381b2SBarry Smith PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, PetscCtx ctx)
101d71ae5a4SJacob Faibussowitsch {
102c4762a1bSJed Brown AppCtx *actx = (AppCtx *)ctx;
103c4762a1bSJed Brown
104c4762a1bSJed Brown PetscFunctionBegin;
1059566063dSJacob Faibussowitsch /* PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); */
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 }
1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
115c4762a1bSJed Brown }
116c4762a1bSJed Brown
117c4762a1bSJed Brown /*
118c4762a1bSJed Brown Defines the ODE passed to the ODE solver
119c4762a1bSJed Brown */
IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,PetscCtx ctx)120*2a8381b2SBarry Smith static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, PetscCtx ctx)
121d71ae5a4SJacob 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));
1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
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 */
IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,PetscCtx ctx)149*2a8381b2SBarry Smith static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
150d71ae5a4SJacob 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 }
1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
183c4762a1bSJed Brown }
184c4762a1bSJed Brown
185c4762a1bSJed Brown /* Matrix JacobianP is constant so that it only needs to be evaluated once */
RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,PetscCtx ctx)186*2a8381b2SBarry Smith static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, PetscCtx ctx)
187d71ae5a4SJacob Faibussowitsch {
188c4762a1bSJed Brown PetscFunctionBeginUser;
1893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
190c4762a1bSJed Brown }
191c4762a1bSJed Brown
main(int argc,char ** argv)192d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
193d71ae5a4SJacob Faibussowitsch {
194c4762a1bSJed Brown TS ts; /* ODE integrator */
195c4762a1bSJed Brown Vec U; /* solution will be stored here */
196c4762a1bSJed Brown Mat A; /* Jacobian matrix */
197c4762a1bSJed Brown Mat Ap; /* dfdp */
198c4762a1bSJed Brown PetscMPIInt size;
199c4762a1bSJed Brown PetscInt n = 2;
200c4762a1bSJed Brown PetscScalar *u, *v;
201c4762a1bSJed Brown AppCtx app;
202c4762a1bSJed Brown PetscInt direction[1];
203c4762a1bSJed Brown PetscBool terminate[1];
204c4762a1bSJed Brown Vec lambda[2], mu[2];
205c4762a1bSJed Brown PetscReal tend;
206c4762a1bSJed Brown
207c4762a1bSJed Brown FILE *f;
208c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209c4762a1bSJed Brown Initialize program
210c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211327415f7SBarry Smith PetscFunctionBeginUser;
212c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2143c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
215c4762a1bSJed Brown app.mode = 1;
216c4762a1bSJed Brown app.lambda1 = 2.75;
217c4762a1bSJed Brown app.lambda2 = 0.36;
218c4762a1bSJed Brown tend = 0.125;
219d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1adj options", "");
220c4762a1bSJed Brown {
2219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL));
2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL));
2239566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tend", "", "", tend, &tend, NULL));
224c4762a1bSJed Brown }
225d0609cedSBarry Smith PetscOptionsEnd();
226c4762a1bSJed Brown
227c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228c4762a1bSJed Brown Create necessary matrix and vectors
229c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2309566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2319566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
2329566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE));
2339566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
2349566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
235c4762a1bSJed Brown
2369566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL));
237c4762a1bSJed Brown
2389566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Ap));
2399566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Ap, n, 1, PETSC_DETERMINE, PETSC_DETERMINE));
2409566063dSJacob Faibussowitsch PetscCall(MatSetType(Ap, MATDENSE));
2419566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Ap));
2429566063dSJacob Faibussowitsch PetscCall(MatSetUp(Ap));
2439566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Ap)); /* initialize to zeros */
244c4762a1bSJed Brown
2459566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u));
246c4762a1bSJed Brown u[0] = 0;
247c4762a1bSJed Brown u[1] = 1;
2489566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u));
249c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
250c4762a1bSJed Brown Create timestepping solver context
251c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2529566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2539566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
2549566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN));
2558434afd1SBarry Smith PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, &app));
2568434afd1SBarry Smith PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobianFn *)IJacobian, &app));
2579566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ts, Ap, RHSJacobianP, &app));
258c4762a1bSJed Brown
259c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
260c4762a1bSJed Brown Set initial conditions
261c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2629566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U));
263c4762a1bSJed Brown
264c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
265c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used
266c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2679566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts));
268c4762a1bSJed Brown
269c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
270c4762a1bSJed Brown Set solver options
271c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2729566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, tend));
2739566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
2749566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 1. / 256.));
2759566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
276c4762a1bSJed Brown
277c4762a1bSJed Brown /* Set directions and terminate flags for the two events */
278c4762a1bSJed Brown direction[0] = 0;
279c4762a1bSJed Brown terminate[0] = PETSC_FALSE;
2809566063dSJacob Faibussowitsch PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app));
281c4762a1bSJed Brown
282c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
283c4762a1bSJed Brown Run timestepping solver
284c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2859566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U));
286c4762a1bSJed Brown
287c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
288c4762a1bSJed Brown Adjoint model starts here
289c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2909566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lambda[0], NULL));
2919566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lambda[1], NULL));
292c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */
2939566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(lambda[0]));
2949566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(lambda[1]));
2959566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &u));
296c4762a1bSJed Brown u[0] = 1.;
2979566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &u));
2989566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[1], &u));
299c4762a1bSJed Brown u[1] = 1.;
3009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[1], &u));
301c4762a1bSJed Brown
3029566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Ap, &mu[0], NULL));
3039566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Ap, &mu[1], NULL));
3049566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mu[0]));
3059566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mu[1]));
3069566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 2, lambda, mu));
307c4762a1bSJed Brown
3089566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts));
309c4762a1bSJed Brown
310c4762a1bSJed Brown /*
3119566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD));
3129566063dSJacob Faibussowitsch PetscCall(VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD));
3139566063dSJacob Faibussowitsch PetscCall(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD));
3149566063dSJacob Faibussowitsch PetscCall(VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD));
315c4762a1bSJed Brown */
3169566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &u));
3179566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[1], &v));
318c4762a1bSJed Brown f = fopen("adj_mu.out", "a");
31963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, f, "%20.15lf %20.15lf %20.15lf\n", (double)tend, (double)PetscRealPart(u[0]), (double)PetscRealPart(v[0])));
3209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &u));
3219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[1], &v));
322c4762a1bSJed Brown fclose(f);
323c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
324c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed.
325c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
3279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U));
3289566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
329c4762a1bSJed Brown
3309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ap));
3319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0]));
3329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[1]));
3339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[0]));
3349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[1]));
3359566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
336b122ec5aSJacob Faibussowitsch return 0;
337c4762a1bSJed Brown }
338c4762a1bSJed Brown
339c4762a1bSJed Brown /*TEST
340c4762a1bSJed Brown
341c4762a1bSJed Brown build:
342c4762a1bSJed Brown requires: !complex
343c4762a1bSJed Brown
344c4762a1bSJed Brown test:
345fe4ad979SIlya Fursov args: -ts_monitor -ts_adjoint_monitor
346c4762a1bSJed Brown
347c4762a1bSJed Brown TEST*/
348