xref: /petsc/src/ts/tutorials/hybrid/ex1adj.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
27c4762a1bSJed Brown PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx)
28c4762a1bSJed Brown {
29c4762a1bSJed Brown   AppCtx            *actx=(AppCtx*)ctx;
30c4762a1bSJed Brown   const PetscScalar *u;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   PetscFunctionBegin;
339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
34c4762a1bSJed Brown   if (actx->mode == 1) {
35c4762a1bSJed Brown     fvalue[0] = u[1]-actx->lambda1*u[0];
36c4762a1bSJed Brown   }else if (actx->mode == 2) {
37c4762a1bSJed Brown     fvalue[0] = u[1]-actx->lambda2*u[0];
38c4762a1bSJed Brown   }
399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
40c4762a1bSJed Brown   PetscFunctionReturn(0);
41c4762a1bSJed Brown }
42c4762a1bSJed Brown 
43c4762a1bSJed Brown PetscErrorCode ShiftGradients(TS ts,Vec U,AppCtx *actx)
44c4762a1bSJed Brown {
45c4762a1bSJed Brown   Vec               *lambda,*mu;
46c4762a1bSJed Brown   PetscScalar       *x,*y;
47c4762a1bSJed Brown   const PetscScalar *u;
48c4762a1bSJed Brown   PetscScalar       tmp[2],A1[2][2],A2[2],denorm;
49c4762a1bSJed Brown   PetscInt          numcost;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   PetscFunctionBegin;
529566063dSJacob Faibussowitsch   PetscCall(TSGetCostGradients(ts,&numcost,&lambda,&mu));
539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   if (actx->mode==2) {
56c4762a1bSJed Brown     denorm = -actx->lambda1*(u[0]-100.*u[1])+1.*(10.*u[0]+u[1]);
57c4762a1bSJed Brown     A1[0][0] = 110.*u[1]*(-actx->lambda1)/denorm+1.;
58c4762a1bSJed Brown     A1[0][1] = -110.*u[0]*(-actx->lambda1)/denorm;
59c4762a1bSJed Brown     A1[1][0] = 110.*u[1]*1./denorm;
60c4762a1bSJed Brown     A1[1][1] = -110.*u[0]*1./denorm+1.;
61c4762a1bSJed Brown 
62c4762a1bSJed Brown     A2[0] = 110.*u[1]*(-u[0])/denorm;
63c4762a1bSJed Brown     A2[1] = -110.*u[0]*(-u[0])/denorm;
64c4762a1bSJed Brown   } else {
65c4762a1bSJed Brown     denorm = -actx->lambda2*(u[0]+10.*u[1])+1.*(-100.*u[0]+u[1]);
66c4762a1bSJed Brown     A1[0][0] = 110.*u[1]*(actx->lambda2)/denorm+1;
67c4762a1bSJed Brown     A1[0][1] = -110.*u[0]*(actx->lambda2)/denorm;
68c4762a1bSJed Brown     A1[1][0] = -110.*u[1]*1./denorm;
69c4762a1bSJed Brown     A1[1][1] = 110.*u[0]*1./denorm+1.;
70c4762a1bSJed Brown 
71c4762a1bSJed Brown     A2[0] = 0;
72c4762a1bSJed Brown     A2[1] = 0;
73c4762a1bSJed Brown   }
74c4762a1bSJed Brown 
759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0],&x));
789566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0],&y));
79c4762a1bSJed Brown   tmp[0] = A1[0][0]*x[0]+A1[0][1]*x[1];
80c4762a1bSJed Brown   tmp[1] = A1[1][0]*x[0]+A1[1][1]*x[1];
81c4762a1bSJed Brown   y[0]   = y[0] + A2[0]*x[0]+A2[1]*x[1];
82c4762a1bSJed Brown   x[0]   = tmp[0];
83c4762a1bSJed Brown   x[1]   = tmp[1];
849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0],&y));
859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0],&x));
86c4762a1bSJed Brown 
879566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[1],&x));
889566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[1],&y));
89c4762a1bSJed Brown   tmp[0] = A1[0][0]*x[0]+A1[0][1]*x[1];
90c4762a1bSJed Brown   tmp[1] = A1[1][0]*x[0]+A1[1][1]*x[1];
91c4762a1bSJed Brown   y[0]   = y[0] + A2[0]*x[0]+A2[1]*x[1];
92c4762a1bSJed Brown   x[0]   = tmp[0];
93c4762a1bSJed Brown   x[1]   = tmp[1];
949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[1],&y));
959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[1],&x));
96c4762a1bSJed Brown   PetscFunctionReturn(0);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown 
99c4762a1bSJed Brown PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx)
100c4762a1bSJed Brown {
101c4762a1bSJed Brown   AppCtx         *actx=(AppCtx*)ctx;
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   PetscFunctionBegin;
1049566063dSJacob Faibussowitsch   /* PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); */
105c4762a1bSJed Brown   if (!forwardsolve) {
1069566063dSJacob Faibussowitsch     PetscCall(ShiftGradients(ts,U,actx));
107c4762a1bSJed Brown   }
108c4762a1bSJed Brown   if (actx->mode == 1) {
109c4762a1bSJed Brown     actx->mode = 2;
1109566063dSJacob Faibussowitsch     /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t)); */
111c4762a1bSJed Brown   } else if (actx->mode == 2) {
112c4762a1bSJed Brown     actx->mode = 1;
1139566063dSJacob Faibussowitsch     /* PetscCall(PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t)); */
114c4762a1bSJed Brown   }
115c4762a1bSJed Brown   PetscFunctionReturn(0);
116c4762a1bSJed Brown }
117c4762a1bSJed Brown 
118c4762a1bSJed Brown /*
119c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
120c4762a1bSJed Brown */
121c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
122c4762a1bSJed Brown {
123c4762a1bSJed Brown   AppCtx            *actx=(AppCtx*)ctx;
124c4762a1bSJed Brown   PetscScalar       *f;
125c4762a1bSJed Brown   const PetscScalar *u,*udot;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   PetscFunctionBegin;
128c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
1319566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   if (actx->mode == 1) {
134c4762a1bSJed Brown     f[0] = udot[0]-u[0]+100*u[1];
135c4762a1bSJed Brown     f[1] = udot[1]-10*u[0]-u[1];
136c4762a1bSJed Brown   } else if (actx->mode == 2) {
137c4762a1bSJed Brown     f[0] = udot[0]-u[0]-10*u[1];
138c4762a1bSJed Brown     f[1] = udot[1]+100*u[0]-u[1];
139c4762a1bSJed Brown   }
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot,&udot));
1439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
144c4762a1bSJed Brown   PetscFunctionReturn(0);
145c4762a1bSJed Brown }
146c4762a1bSJed Brown 
147c4762a1bSJed Brown /*
148c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
149c4762a1bSJed Brown */
150c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx)
151c4762a1bSJed Brown {
152c4762a1bSJed Brown   AppCtx            *actx=(AppCtx*)ctx;
153c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
154c4762a1bSJed Brown   PetscScalar       J[2][2];
155c4762a1bSJed Brown   const PetscScalar *u,*udot;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   PetscFunctionBegin;
1589566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1599566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   if (actx->mode == 1) {
162c4762a1bSJed Brown     J[0][0] = a-1;                       J[0][1] = 100;
163c4762a1bSJed Brown     J[1][0] = -10;                       J[1][1] = a-1;
164c4762a1bSJed Brown   } else if (actx->mode == 2) {
165c4762a1bSJed Brown     J[0][0] = a-1;                       J[0][1] = -10;
166c4762a1bSJed Brown     J[1][0] = 100;                       J[1][1] = a-1;
167c4762a1bSJed Brown   }
1689566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
169c4762a1bSJed Brown 
1709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot,&udot));
172c4762a1bSJed Brown 
1739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
175c4762a1bSJed Brown   if (A != B) {
1769566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1779566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
178c4762a1bSJed Brown   }
179c4762a1bSJed Brown   PetscFunctionReturn(0);
180c4762a1bSJed Brown }
181c4762a1bSJed Brown 
182c4762a1bSJed Brown /* Matrix JacobianP is constant so that it only needs to be evaluated once */
183c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A, void *ctx)
184c4762a1bSJed Brown {
185c4762a1bSJed Brown   PetscFunctionBeginUser;
186c4762a1bSJed Brown   PetscFunctionReturn(0);
187c4762a1bSJed Brown }
188c4762a1bSJed Brown 
189c4762a1bSJed Brown int main(int argc,char **argv)
190c4762a1bSJed Brown {
191c4762a1bSJed Brown   TS             ts;            /* ODE integrator */
192c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
193c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
194c4762a1bSJed Brown   Mat            Ap;            /* dfdp */
195c4762a1bSJed Brown   PetscMPIInt    size;
196c4762a1bSJed Brown   PetscInt       n = 2;
197c4762a1bSJed Brown   PetscScalar    *u,*v;
198c4762a1bSJed Brown   AppCtx         app;
199c4762a1bSJed Brown   PetscInt       direction[1];
200c4762a1bSJed Brown   PetscBool      terminate[1];
201c4762a1bSJed Brown   Vec            lambda[2],mu[2];
202c4762a1bSJed Brown   PetscReal      tend;
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   FILE           *f;
205c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206c4762a1bSJed Brown      Initialize program
207c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
208*327415f7SBarry Smith   PetscFunctionBeginUser;
2099566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
2109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
2113c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
212c4762a1bSJed Brown   app.mode = 1;
213c4762a1bSJed Brown   app.lambda1 = 2.75;
214c4762a1bSJed Brown   app.lambda2 = 0.36;
215c4762a1bSJed Brown   tend = 0.125;
216d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex1adj options","");
217c4762a1bSJed Brown   {
2189566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda1","","",app.lambda1,&app.lambda1,NULL));
2199566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-lambda2","","",app.lambda2,&app.lambda2,NULL));
2209566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tend","","",tend,&tend,NULL));
221c4762a1bSJed Brown   }
222d0609cedSBarry Smith   PetscOptionsEnd();
223c4762a1bSJed Brown 
224c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225c4762a1bSJed Brown     Create necessary matrix and vectors
226c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2279566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
2289566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
2299566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATDENSE));
2309566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2319566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
232c4762a1bSJed Brown 
2339566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&U,NULL));
234c4762a1bSJed Brown 
2359566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&Ap));
2369566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Ap,n,1,PETSC_DETERMINE,PETSC_DETERMINE));
2379566063dSJacob Faibussowitsch   PetscCall(MatSetType(Ap,MATDENSE));
2389566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Ap));
2399566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Ap));
2409566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(Ap)); /* initialize to zeros */
241c4762a1bSJed Brown 
2429566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U,&u));
243c4762a1bSJed Brown   u[0] = 0;
244c4762a1bSJed Brown   u[1] = 1;
2459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U,&u));
246c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
247c4762a1bSJed Brown      Create timestepping solver context
248c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2499566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
2509566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
2519566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSCN));
2529566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,(TSIFunction)IFunction,&app));
2539566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&app));
2549566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ts,Ap,RHSJacobianP,&app));
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
257c4762a1bSJed Brown      Set initial conditions
258c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2599566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,U));
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
263c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2649566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267c4762a1bSJed Brown      Set solver options
268c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2699566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,tend));
2709566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
2719566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,1./256.));
2729566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
273c4762a1bSJed Brown 
274c4762a1bSJed Brown   /* Set directions and terminate flags for the two events */
275c4762a1bSJed Brown   direction[0] = 0;
276c4762a1bSJed Brown   terminate[0] = PETSC_FALSE;
2779566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts,1,direction,terminate,EventFunction,PostEventFunction,(void*)&app));
278c4762a1bSJed Brown 
279c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
280c4762a1bSJed Brown      Run timestepping solver
281c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2829566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,U));
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
285c4762a1bSJed Brown      Adjoint model starts here
286c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2879566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&lambda[0],NULL));
2889566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&lambda[1],NULL));
289c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
2909566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(lambda[0]));
2919566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(lambda[1]));
2929566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0],&u));
293c4762a1bSJed Brown   u[0] = 1.;
2949566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0],&u));
2959566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[1],&u));
296c4762a1bSJed Brown   u[1] = 1.;
2979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[1],&u));
298c4762a1bSJed Brown 
2999566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Ap,&mu[0],NULL));
3009566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Ap,&mu[1],NULL));
3019566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(mu[0]));
3029566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(mu[1]));
3039566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts,2,lambda,mu));
304c4762a1bSJed Brown 
3059566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
306c4762a1bSJed Brown 
307c4762a1bSJed Brown   /*
3089566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD));
3099566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD));
3109566063dSJacob Faibussowitsch   PetscCall(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD));
3119566063dSJacob Faibussowitsch   PetscCall(VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD));
312c4762a1bSJed Brown   */
3139566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0],&u));
3149566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[1],&v));
315c4762a1bSJed Brown   f = fopen("adj_mu.out", "a");
31663a3b9bcSJacob Faibussowitsch   PetscCall(PetscFPrintf(PETSC_COMM_WORLD,f,"%20.15lf %20.15lf %20.15lf\n",(double)tend,(double)PetscRealPart(u[0]),(double)PetscRealPart(v[0])));
3179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0],&u));
3189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[1],&v));
319c4762a1bSJed Brown   fclose(f);
320c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
321c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
322c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3249566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3259566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
326c4762a1bSJed Brown 
3279566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ap));
3289566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[0]));
3299566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[1]));
3309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[0]));
3319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[1]));
3329566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
333b122ec5aSJacob Faibussowitsch   return 0;
334c4762a1bSJed Brown }
335c4762a1bSJed Brown 
336c4762a1bSJed Brown /*TEST
337c4762a1bSJed Brown 
338c4762a1bSJed Brown    build:
339c4762a1bSJed Brown       requires: !complex
340c4762a1bSJed Brown 
341c4762a1bSJed Brown    test:
342c4762a1bSJed Brown       args: -ts_monitor -ts_adjoint_monitor
343c4762a1bSJed Brown 
344c4762a1bSJed Brown TEST*/
345