xref: /petsc/src/ts/interface/sensitivity/tssen.c (revision 5bd1e5768a3b141382a3ebb2e780cd2d3b3bfbf0)
1 #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2 #include <petscdraw.h>
3 
4 PetscLogEvent TS_AdjointStep,TS_ForwardStep,TS_JacobianPEval;
5 
6 /* #define TSADJOINT_STAGE */
7 
8 /* ------------------------ Sensitivity Context ---------------------------*/
9 
10 /*@C
11   TSSetRHSJacobianP - Sets the function that computes the Jacobian of G w.r.t. the parameters P where U_t = G(U,P,t), as well as the location to store the matrix.
12 
13   Logically Collective on TS
14 
15   Input Parameters:
16 + ts - TS context obtained from TSCreate()
17 . Amat - JacobianP matrix
18 . func - function
19 - ctx - [optional] user-defined function context
20 
21   Calling sequence of func:
22 $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
23 +   t - current timestep
24 .   U - input vector (current ODE solution)
25 .   A - output matrix
26 -   ctx - [optional] user-defined function context
27 
28   Level: intermediate
29 
30   Notes:
31     Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p
32 
33 .seealso: TSGetRHSJacobianP()
34 @*/
35 PetscErrorCode TSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
36 {
37   PetscErrorCode ierr;
38 
39   PetscFunctionBegin;
40   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
41   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
42 
43   ts->rhsjacobianp    = func;
44   ts->rhsjacobianpctx = ctx;
45   if(Amat) {
46     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
47     ierr = MatDestroy(&ts->Jacprhs);CHKERRQ(ierr);
48     ts->Jacprhs = Amat;
49   }
50   PetscFunctionReturn(0);
51 }
52 
53 /*@C
54   TSGetRHSJacobianP - Gets the function that computes the Jacobian of G w.r.t. the parameters P where U_t = G(U,P,t), as well as the location to store the matrix.
55 
56   Logically Collective on TS
57 
58   Input Parameters:
59 . ts - TS context obtained from TSCreate()
60 
61   Output Parameters:
62 + Amat - JacobianP matrix
63 . func - function
64 - ctx - [optional] user-defined function context
65 
66   Calling sequence of func:
67 $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
68 +   t - current timestep
69 .   U - input vector (current ODE solution)
70 .   A - output matrix
71 -   ctx - [optional] user-defined function context
72 
73   Level: intermediate
74 
75   Notes:
76     Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p
77 
78 .seealso: TSSetRHSJacobianP()
79 @*/
80 PetscErrorCode TSGetRHSJacobianP(TS ts,Mat *Amat,PetscErrorCode (**func)(TS,PetscReal,Vec,Mat,void*),void **ctx)
81 {
82   PetscFunctionBegin;
83   if (func) *func = ts->rhsjacobianp;
84   if (ctx) *ctx  = ts->rhsjacobianpctx;
85   if (Amat) *Amat = ts->Jacprhs;
86   PetscFunctionReturn(0);
87 }
88 
89 /*@C
90   TSComputeRHSJacobianP - Runs the user-defined JacobianP function.
91 
92   Collective on TS
93 
94   Input Parameters:
95 . ts   - The TS context obtained from TSCreate()
96 
97   Level: developer
98 
99 .seealso: TSSetRHSJacobianP()
100 @*/
101 PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec U,Mat Amat)
102 {
103   PetscErrorCode ierr;
104 
105   PetscFunctionBegin;
106   if (!Amat) PetscFunctionReturn(0);
107   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
108   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
109 
110   PetscStackPush("TS user JacobianP function for sensitivity analysis");
111   ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx);CHKERRQ(ierr);
112   PetscStackPop;
113   PetscFunctionReturn(0);
114 }
115 
116 /*@C
117   TSSetIJacobianP - Sets the function that computes the Jacobian of F w.r.t. the parameters P where F(Udot,U,t) = G(U,P,t), as well as the location to store the matrix.
118 
119   Logically Collective on TS
120 
121   Input Parameters:
122 + ts - TS context obtained from TSCreate()
123 . Amat - JacobianP matrix
124 . func - function
125 - ctx - [optional] user-defined function context
126 
127   Calling sequence of func:
128 $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
129 +   t - current timestep
130 .   U - input vector (current ODE solution)
131 .   Udot - time derivative of state vector
132 .   shift - shift to apply, see note below
133 .   A - output matrix
134 -   ctx - [optional] user-defined function context
135 
136   Level: intermediate
137 
138   Notes:
139     Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p
140 
141 .seealso:
142 @*/
143 PetscErrorCode TSSetIJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void *ctx)
144 {
145   PetscErrorCode ierr;
146 
147   PetscFunctionBegin;
148   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
149   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
150 
151   ts->ijacobianp    = func;
152   ts->ijacobianpctx = ctx;
153   if(Amat) {
154     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
155     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
156     ts->Jacp = Amat;
157   }
158   PetscFunctionReturn(0);
159 }
160 
161 /*@C
162   TSComputeIJacobianP - Runs the user-defined IJacobianP function.
163 
164   Collective on TS
165 
166   Input Parameters:
167 + ts - the TS context
168 . t - current timestep
169 . U - state vector
170 . Udot - time derivative of state vector
171 . shift - shift to apply, see note below
172 - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate
173 
174   Output Parameters:
175 . A - Jacobian matrix
176 
177   Level: developer
178 
179 .seealso: TSSetIJacobianP()
180 @*/
181 PetscErrorCode TSComputeIJacobianP(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat Amat,PetscBool imex)
182 {
183   PetscErrorCode ierr;
184 
185   PetscFunctionBegin;
186   if (!Amat) PetscFunctionReturn(0);
187   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
188   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
189   PetscValidHeaderSpecific(Udot,VEC_CLASSID,4);
190 
191   ierr = PetscLogEventBegin(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr);
192   if (ts->ijacobianp) {
193     PetscStackPush("TS user JacobianP function for sensitivity analysis");
194     ierr = (*ts->ijacobianp)(ts,t,U,Udot,shift,Amat,ts->ijacobianpctx);CHKERRQ(ierr);
195     PetscStackPop;
196   }
197   if (imex) {
198     if (!ts->ijacobianp) {  /* system was written as Udot = G(t,U) */
199       PetscBool assembled;
200       ierr = MatZeroEntries(Amat);CHKERRQ(ierr);
201       ierr = MatAssembled(Amat,&assembled);CHKERRQ(ierr);
202       if (!assembled) {
203         ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
204         ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
205       }
206     }
207   } else {
208     if (ts->rhsjacobianp) {
209       ierr = TSComputeRHSJacobianP(ts,t,U,ts->Jacprhs);CHKERRQ(ierr);
210     }
211     if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */
212       ierr = MatScale(Amat,-1);CHKERRQ(ierr);
213     } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */
214       MatStructure axpy = DIFFERENT_NONZERO_PATTERN;
215       if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */
216         ierr = MatZeroEntries(Amat);CHKERRQ(ierr);
217       }
218       ierr = MatAXPY(Amat,-1,ts->Jacprhs,axpy);CHKERRQ(ierr);
219     }
220   }
221   ierr = PetscLogEventEnd(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr);
222   PetscFunctionReturn(0);
223 }
224 
225 /*@C
226     TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions
227 
228     Logically Collective on TS
229 
230     Input Parameters:
231 +   ts - the TS context obtained from TSCreate()
232 .   numcost - number of gradients to be computed, this is the number of cost functions
233 .   costintegral - vector that stores the integral values
234 .   rf - routine for evaluating the integrand function
235 .   drduf - function that computes the gradients of the r's with respect to u
236 .   drdpf - function that computes the gradients of the r's with respect to p, can be NULL if parametric sensitivity is not desired (mu=NULL)
237 .   fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
238 -   ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
239 
240     Calling sequence of rf:
241 $   PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx);
242 
243     Calling sequence of drduf:
244 $   PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx);
245 
246     Calling sequence of drdpf:
247 $   PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx);
248 
249     Level: deprecated
250 
251     Notes:
252     For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions
253 
254 .seealso: TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients()
255 @*/
256 PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*),
257                                                           PetscErrorCode (*drduf)(TS,PetscReal,Vec,Vec*,void*),
258                                                           PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*),
259                                                           PetscBool fwd,void *ctx)
260 {
261   PetscErrorCode ierr;
262 
263   PetscFunctionBegin;
264   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
265   if (costintegral) PetscValidHeaderSpecific(costintegral,VEC_CLASSID,3);
266   if (ts->numcost && ts->numcost!=numcost) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostGradients() or TSForwardSetIntegralGradients()");
267   if (!ts->numcost) ts->numcost=numcost;
268 
269   if (costintegral) {
270     ierr = PetscObjectReference((PetscObject)costintegral);CHKERRQ(ierr);
271     ierr = VecDestroy(&ts->vec_costintegral);CHKERRQ(ierr);
272     ts->vec_costintegral = costintegral;
273   } else {
274     if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */
275       ierr = VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);CHKERRQ(ierr);
276     } else {
277       ierr = VecSet(ts->vec_costintegral,0.0);CHKERRQ(ierr);
278     }
279   }
280   if (!ts->vec_costintegrand) {
281     ierr = VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);CHKERRQ(ierr);
282   } else {
283     ierr = VecSet(ts->vec_costintegrand,0.0);CHKERRQ(ierr);
284   }
285   ts->costintegralfwd  = fwd; /* Evaluate the cost integral in forward run if fwd is true */
286   ts->costintegrand    = rf;
287   ts->costintegrandctx = ctx;
288   ts->drdufunction     = drduf;
289   ts->drdpfunction     = drdpf;
290   PetscFunctionReturn(0);
291 }
292 
293 /*@C
294    TSGetCostIntegral - Returns the values of the integral term in the cost functions.
295    It is valid to call the routine after a backward run.
296 
297    Not Collective
298 
299    Input Parameter:
300 .  ts - the TS context obtained from TSCreate()
301 
302    Output Parameter:
303 .  v - the vector containing the integrals for each cost function
304 
305    Level: intermediate
306 
307 .seealso: TSSetCostIntegrand()
308 
309 @*/
310 PetscErrorCode  TSGetCostIntegral(TS ts,Vec *v)
311 {
312   TS             quadts;
313   PetscErrorCode ierr;
314 
315   PetscFunctionBegin;
316   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
317   PetscValidPointer(v,2);
318   ierr = TSGetQuadratureTS(ts,NULL,&quadts);CHKERRQ(ierr);
319   *v = quadts->vec_sol;
320   PetscFunctionReturn(0);
321 }
322 
323 /*@C
324    TSComputeCostIntegrand - Evaluates the integral function in the cost functions.
325 
326    Input Parameters:
327 +  ts - the TS context
328 .  t - current time
329 -  U - state vector, i.e. current solution
330 
331    Output Parameter:
332 .  Q - vector of size numcost to hold the outputs
333 
334    Notes:
335    Most users should not need to explicitly call this routine, as it
336    is used internally within the sensitivity analysis context.
337 
338    Level: deprecated
339 
340 .seealso: TSSetCostIntegrand()
341 @*/
342 PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec U,Vec Q)
343 {
344   PetscErrorCode ierr;
345 
346   PetscFunctionBegin;
347   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
348   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
349   PetscValidHeaderSpecific(Q,VEC_CLASSID,4);
350 
351   ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
352   if (ts->costintegrand) {
353     PetscStackPush("TS user integrand in the cost function");
354     ierr = (*ts->costintegrand)(ts,t,U,Q,ts->costintegrandctx);CHKERRQ(ierr);
355     PetscStackPop;
356   } else {
357     ierr = VecZeroEntries(Q);CHKERRQ(ierr);
358   }
359 
360   ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
361   PetscFunctionReturn(0);
362 }
363 
364 /*@C
365   TSComputeDRDUFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian()
366 
367   Level: deprecated
368 
369 @*/
370 PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
371 {
372   PetscErrorCode ierr;
373 
374   PetscFunctionBegin;
375   if (!DRDU) PetscFunctionReturn(0);
376   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
377   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
378 
379   PetscStackPush("TS user DRDU function for sensitivity analysis");
380   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
381   PetscStackPop;
382   PetscFunctionReturn(0);
383 }
384 
385 /*@C
386   TSComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP()
387 
388   Level: deprecated
389 
390 @*/
391 PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
392 {
393   PetscErrorCode ierr;
394 
395   PetscFunctionBegin;
396   if (!DRDP) PetscFunctionReturn(0);
397   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
398   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
399 
400   PetscStackPush("TS user DRDP function for sensitivity analysis");
401   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
402   PetscStackPop;
403   PetscFunctionReturn(0);
404 }
405 
406 /*@C
407   TSSetIHessianProduct - Sets the function that computes the vector-Hessian-vector product. The Hessian is the second-order derivative of F (IFunction) w.r.t. the state variable.
408 
409   Logically Collective on TS
410 
411   Input Parameters:
412 + ts - TS context obtained from TSCreate()
413 . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU
414 . hessianproductfunc1 - vector-Hessian-vector product function for F_UU
415 . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP
416 . hessianproductfunc2 - vector-Hessian-vector product function for F_UP
417 . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU
418 . hessianproductfunc3 - vector-Hessian-vector product function for F_PU
419 . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP
420 . hessianproductfunc4 - vector-Hessian-vector product function for F_PP
421 
422   Calling sequence of ihessianproductfunc:
423 $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx);
424 +   t - current timestep
425 .   U - input vector (current ODE solution)
426 .   Vl - an array of input vectors to be left-multiplied with the Hessian
427 .   Vr - input vector to be right-multiplied with the Hessian
428 .   VHV - an array of output vectors for vector-Hessian-vector product
429 -   ctx - [optional] user-defined function context
430 
431   Level: intermediate
432 
433   Notes:
434   The first Hessian function and the working array are required.
435   As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
436   $ Vl_n^T*F_UP*Vr
437   where the vector Vl_n (n-th element in the array Vl) and Vr are of size N and M respectively, and the Hessian F_UP is of size N x N x M.
438   Each entry of F_UP corresponds to the derivative
439   $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}.
440   The result of the vector-Hessian-vector product for Vl_n needs to be stored in vector VHV_n with the j-th entry being
441   $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]}
442   If the cost function is a scalar, there will be only one vector in Vl and VHV.
443 
444 .seealso:
445 @*/
446 PetscErrorCode TSSetIHessianProduct(TS ts,Vec *ihp1,PetscErrorCode (*ihessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
447                                           Vec *ihp2,PetscErrorCode (*ihessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
448                                           Vec *ihp3,PetscErrorCode (*ihessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
449                                           Vec *ihp4,PetscErrorCode (*ihessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
450                                     void *ctx)
451 {
452   PetscFunctionBegin;
453   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
454   PetscValidPointer(ihp1,2);
455 
456   ts->ihessianproductctx = ctx;
457   if (ihp1) ts->vecs_fuu = ihp1;
458   if (ihp2) ts->vecs_fup = ihp2;
459   if (ihp3) ts->vecs_fpu = ihp3;
460   if (ihp4) ts->vecs_fpp = ihp4;
461   ts->ihessianproduct_fuu = ihessianproductfunc1;
462   ts->ihessianproduct_fup = ihessianproductfunc2;
463   ts->ihessianproduct_fpu = ihessianproductfunc3;
464   ts->ihessianproduct_fpp = ihessianproductfunc4;
465   PetscFunctionReturn(0);
466 }
467 
468 /*@C
469   TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu.
470 
471   Collective on TS
472 
473   Input Parameters:
474 . ts   - The TS context obtained from TSCreate()
475 
476   Notes:
477   TSComputeIHessianProductFunctionUU() is typically used for sensitivity implementation,
478   so most users would not generally call this routine themselves.
479 
480   Level: developer
481 
482 .seealso: TSSetIHessianProduct()
483 @*/
484 PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
485 {
486   PetscErrorCode ierr;
487 
488   PetscFunctionBegin;
489   if (!VHV) PetscFunctionReturn(0);
490   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
491   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
492 
493   if (ts->ihessianproduct_fuu) {
494     PetscStackPush("TS user IHessianProduct function 1 for sensitivity analysis");
495     ierr = (*ts->ihessianproduct_fuu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
496     PetscStackPop;
497   }
498   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
499   if (ts->rhshessianproduct_guu) {
500     PetscInt nadj;
501     ierr = TSComputeRHSHessianProductFunctionUU(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
502     for (nadj=0; nadj<ts->numcost; nadj++) {
503       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
504     }
505   }
506   PetscFunctionReturn(0);
507 }
508 
509 /*@C
510   TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup.
511 
512   Collective on TS
513 
514   Input Parameters:
515 . ts   - The TS context obtained from TSCreate()
516 
517   Notes:
518   TSComputeIHessianProductFunctionUP() is typically used for sensitivity implementation,
519   so most users would not generally call this routine themselves.
520 
521   Level: developer
522 
523 .seealso: TSSetIHessianProduct()
524 @*/
525 PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
526 {
527   PetscErrorCode ierr;
528 
529   PetscFunctionBegin;
530   if (!VHV) PetscFunctionReturn(0);
531   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
532   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
533 
534   if (ts->ihessianproduct_fup) {
535     PetscStackPush("TS user IHessianProduct function 2 for sensitivity analysis");
536     ierr = (*ts->ihessianproduct_fup)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
537     PetscStackPop;
538   }
539   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
540   if (ts->rhshessianproduct_gup) {
541     PetscInt nadj;
542     ierr = TSComputeRHSHessianProductFunctionUP(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
543     for (nadj=0; nadj<ts->numcost; nadj++) {
544       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
545     }
546   }
547   PetscFunctionReturn(0);
548 }
549 
550 /*@C
551   TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu.
552 
553   Collective on TS
554 
555   Input Parameters:
556 . ts   - The TS context obtained from TSCreate()
557 
558   Notes:
559   TSComputeIHessianProductFunctionPU() is typically used for sensitivity implementation,
560   so most users would not generally call this routine themselves.
561 
562   Level: developer
563 
564 .seealso: TSSetIHessianProduct()
565 @*/
566 PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
567 {
568   PetscErrorCode ierr;
569 
570   PetscFunctionBegin;
571   if (!VHV) PetscFunctionReturn(0);
572   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
573   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
574 
575   if (ts->ihessianproduct_fpu) {
576     PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
577     ierr = (*ts->ihessianproduct_fpu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
578     PetscStackPop;
579   }
580   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
581   if (ts->rhshessianproduct_gpu) {
582     PetscInt nadj;
583     ierr = TSComputeRHSHessianProductFunctionPU(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
584     for (nadj=0; nadj<ts->numcost; nadj++) {
585       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
586     }
587   }
588   PetscFunctionReturn(0);
589 }
590 
591 /*@C
592   TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp.
593 
594   Collective on TS
595 
596   Input Parameters:
597 . ts   - The TS context obtained from TSCreate()
598 
599   Notes:
600   TSComputeIHessianProductFunctionPP() is typically used for sensitivity implementation,
601   so most users would not generally call this routine themselves.
602 
603   Level: developer
604 
605 .seealso: TSSetIHessianProduct()
606 @*/
607 PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
608 {
609   PetscErrorCode ierr;
610 
611   PetscFunctionBegin;
612   if (!VHV) PetscFunctionReturn(0);
613   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
614   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
615 
616   if (ts->ihessianproduct_fpp) {
617     PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis");
618     ierr = (*ts->ihessianproduct_fpp)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr);
619     PetscStackPop;
620   }
621   /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */
622   if (ts->rhshessianproduct_gpp) {
623     PetscInt nadj;
624     ierr = TSComputeRHSHessianProductFunctionPP(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr);
625     for (nadj=0; nadj<ts->numcost; nadj++) {
626       ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr);
627     }
628   }
629   PetscFunctionReturn(0);
630 }
631 
632 /*@C
633   TSSetRHSHessianProduct - Sets the function that computes the vector-Hessian-vector product. The Hessian is the second-order derivative of G (RHSFunction) w.r.t. the state variable.
634 
635   Logically Collective on TS
636 
637   Input Parameters:
638 + ts - TS context obtained from TSCreate()
639 . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU
640 . hessianproductfunc1 - vector-Hessian-vector product function for G_UU
641 . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP
642 . hessianproductfunc2 - vector-Hessian-vector product function for G_UP
643 . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU
644 . hessianproductfunc3 - vector-Hessian-vector product function for G_PU
645 . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP
646 . hessianproductfunc4 - vector-Hessian-vector product function for G_PP
647 
648   Calling sequence of ihessianproductfunc:
649 $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx);
650 +   t - current timestep
651 .   U - input vector (current ODE solution)
652 .   Vl - an array of input vectors to be left-multiplied with the Hessian
653 .   Vr - input vector to be right-multiplied with the Hessian
654 .   VHV - an array of output vectors for vector-Hessian-vector product
655 -   ctx - [optional] user-defined function context
656 
657   Level: intermediate
658 
659   Notes:
660   The first Hessian function and the working array are required.
661   As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product
662   $ Vl_n^T*G_UP*Vr
663   where the vector Vl_n (n-th element in the array Vl) and Vr are of size N and M respectively, and the Hessian G_UP is of size N x N x M.
664   Each entry of G_UP corresponds to the derivative
665   $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}.
666   The result of the vector-Hessian-vector product for Vl_n needs to be stored in vector VHV_n with j-th entry being
667   $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]}
668   If the cost function is a scalar, there will be only one vector in Vl and VHV.
669 
670 .seealso:
671 @*/
672 PetscErrorCode TSSetRHSHessianProduct(TS ts,Vec *rhshp1,PetscErrorCode (*rhshessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
673                                           Vec *rhshp2,PetscErrorCode (*rhshessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
674                                           Vec *rhshp3,PetscErrorCode (*rhshessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
675                                           Vec *rhshp4,PetscErrorCode (*rhshessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*),
676                                     void *ctx)
677 {
678   PetscFunctionBegin;
679   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
680   PetscValidPointer(rhshp1,2);
681 
682   ts->rhshessianproductctx = ctx;
683   if (rhshp1) ts->vecs_guu = rhshp1;
684   if (rhshp2) ts->vecs_gup = rhshp2;
685   if (rhshp3) ts->vecs_gpu = rhshp3;
686   if (rhshp4) ts->vecs_gpp = rhshp4;
687   ts->rhshessianproduct_guu = rhshessianproductfunc1;
688   ts->rhshessianproduct_gup = rhshessianproductfunc2;
689   ts->rhshessianproduct_gpu = rhshessianproductfunc3;
690   ts->rhshessianproduct_gpp = rhshessianproductfunc4;
691   PetscFunctionReturn(0);
692 }
693 
694 /*@C
695   TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu.
696 
697   Collective on TS
698 
699   Input Parameters:
700 . ts   - The TS context obtained from TSCreate()
701 
702   Notes:
703   TSComputeRHSHessianProductFunctionUU() is typically used for sensitivity implementation,
704   so most users would not generally call this routine themselves.
705 
706   Level: developer
707 
708 .seealso: TSSetRHSHessianProduct()
709 @*/
710 PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
711 {
712   PetscErrorCode ierr;
713 
714   PetscFunctionBegin;
715   if (!VHV) PetscFunctionReturn(0);
716   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
717   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
718 
719   PetscStackPush("TS user RHSHessianProduct function 1 for sensitivity analysis");
720   ierr = (*ts->rhshessianproduct_guu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
721   PetscStackPop;
722   PetscFunctionReturn(0);
723 }
724 
725 /*@C
726   TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup.
727 
728   Collective on TS
729 
730   Input Parameters:
731 . ts   - The TS context obtained from TSCreate()
732 
733   Notes:
734   TSComputeRHSHessianProductFunctionUP() is typically used for sensitivity implementation,
735   so most users would not generally call this routine themselves.
736 
737   Level: developer
738 
739 .seealso: TSSetRHSHessianProduct()
740 @*/
741 PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
742 {
743   PetscErrorCode ierr;
744 
745   PetscFunctionBegin;
746   if (!VHV) PetscFunctionReturn(0);
747   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
748   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
749 
750   PetscStackPush("TS user RHSHessianProduct function 2 for sensitivity analysis");
751   ierr = (*ts->rhshessianproduct_gup)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
752   PetscStackPop;
753   PetscFunctionReturn(0);
754 }
755 
756 /*@C
757   TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu.
758 
759   Collective on TS
760 
761   Input Parameters:
762 . ts   - The TS context obtained from TSCreate()
763 
764   Notes:
765   TSComputeRHSHessianProductFunctionPU() is typically used for sensitivity implementation,
766   so most users would not generally call this routine themselves.
767 
768   Level: developer
769 
770 .seealso: TSSetRHSHessianProduct()
771 @*/
772 PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
773 {
774   PetscErrorCode ierr;
775 
776   PetscFunctionBegin;
777   if (!VHV) PetscFunctionReturn(0);
778   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
779   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
780 
781   PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis");
782   ierr = (*ts->rhshessianproduct_gpu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
783   PetscStackPop;
784   PetscFunctionReturn(0);
785 }
786 
787 /*@C
788   TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp.
789 
790   Collective on TS
791 
792   Input Parameters:
793 . ts   - The TS context obtained from TSCreate()
794 
795   Notes:
796   TSComputeRHSHessianProductFunctionPP() is typically used for sensitivity implementation,
797   so most users would not generally call this routine themselves.
798 
799   Level: developer
800 
801 .seealso: TSSetRHSHessianProduct()
802 @*/
803 PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV)
804 {
805   PetscErrorCode ierr;
806 
807   PetscFunctionBegin;
808   if (!VHV) PetscFunctionReturn(0);
809   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
810   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
811 
812   PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis");
813   ierr = (*ts->rhshessianproduct_gpp)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr);
814   PetscStackPop;
815   PetscFunctionReturn(0);
816 }
817 
818 /* --------------------------- Adjoint sensitivity ---------------------------*/
819 
820 /*@
821    TSSetCostGradients - Sets the initial value of the gradients of the cost function w.r.t. initial values and w.r.t. the problem parameters
822       for use by the TSAdjoint routines.
823 
824    Logically Collective on TS
825 
826    Input Parameters:
827 +  ts - the TS context obtained from TSCreate()
828 .  lambda - gradients with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector
829 -  mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
830 
831    Level: beginner
832 
833    Notes:
834     the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime  mu_i = df/dp|finaltime
835 
836    After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities
837 
838 .seealso TSGetCostGradients()
839 @*/
840 PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu)
841 {
842   PetscFunctionBegin;
843   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
844   PetscValidPointer(lambda,2);
845   ts->vecs_sensi  = lambda;
846   ts->vecs_sensip = mu;
847   if (ts->numcost && ts->numcost!=numcost) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand");
848   ts->numcost  = numcost;
849   PetscFunctionReturn(0);
850 }
851 
852 /*@
853    TSGetCostGradients - Returns the gradients from the TSAdjointSolve()
854 
855    Not Collective, but Vec returned is parallel if TS is parallel
856 
857    Input Parameter:
858 .  ts - the TS context obtained from TSCreate()
859 
860    Output Parameter:
861 +  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
862 -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
863 
864    Level: intermediate
865 
866 .seealso: TSSetCostGradients()
867 @*/
868 PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu)
869 {
870   PetscFunctionBegin;
871   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
872   if (numcost) *numcost = ts->numcost;
873   if (lambda)  *lambda  = ts->vecs_sensi;
874   if (mu)      *mu      = ts->vecs_sensip;
875   PetscFunctionReturn(0);
876 }
877 
878 /*@
879    TSSetCostHessianProducts - Sets the initial value of the Hessian-vector products of the cost function w.r.t. initial values and w.r.t. the problem parameters
880       for use by the TSAdjoint routines.
881 
882    Logically Collective on TS
883 
884    Input Parameters:
885 +  ts - the TS context obtained from TSCreate()
886 .  numcost - number of cost functions
887 .  lambda2 - Hessian-vector product with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector
888 .  mu2 - Hessian-vector product with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
889 -  dir - the direction vector that are multiplied with the Hessian of the cost functions
890 
891    Level: beginner
892 
893    Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system
894 
895    For second-order adjoint, one needs to call this function and then TSAdjointSetForward() before TSSolve().
896 
897    After TSAdjointSolve() is called, the lamba2 and the mu2 will contain the computed second-order adjoint sensitivities, and can be used to produce Hessian-vector product (not the full Hessian matrix). Users must provide a direction vector; it is usually generated by an optimization solver.
898 
899    Passing NULL for lambda2 disables the second-order calculation.
900 .seealso: TSAdjointSetForward()
901 @*/
902 PetscErrorCode TSSetCostHessianProducts(TS ts,PetscInt numcost,Vec *lambda2,Vec *mu2,Vec dir)
903 {
904   PetscFunctionBegin;
905   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
906   if (ts->numcost && ts->numcost!=numcost) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand");
907   ts->numcost       = numcost;
908   ts->vecs_sensi2   = lambda2;
909   ts->vecs_sensi2p  = mu2;
910   ts->vec_dir       = dir;
911   PetscFunctionReturn(0);
912 }
913 
914 /*@
915    TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve()
916 
917    Not Collective, but Vec returned is parallel if TS is parallel
918 
919    Input Parameter:
920 .  ts - the TS context obtained from TSCreate()
921 
922    Output Parameter:
923 +  numcost - number of cost functions
924 .  lambda2 - Hessian-vector product with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector
925 .  mu2 - Hessian-vector product with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
926 -  dir - the direction vector that are multiplied with the Hessian of the cost functions
927 
928    Level: intermediate
929 
930 .seealso: TSSetCostHessianProducts()
931 @*/
932 PetscErrorCode TSGetCostHessianProducts(TS ts,PetscInt *numcost,Vec **lambda2,Vec **mu2, Vec *dir)
933 {
934   PetscFunctionBegin;
935   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
936   if (numcost) *numcost = ts->numcost;
937   if (lambda2) *lambda2 = ts->vecs_sensi2;
938   if (mu2)     *mu2     = ts->vecs_sensi2p;
939   if (dir)     *dir     = ts->vec_dir;
940   PetscFunctionReturn(0);
941 }
942 
943 /*@
944   TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities
945 
946   Logically Collective on TS
947 
948   Input Parameters:
949 +  ts - the TS context obtained from TSCreate()
950 -  didp - the derivative of initial values w.r.t. parameters
951 
952   Level: intermediate
953 
954   Notes: When computing sensitivies w.r.t. initial condition, set didp to NULL so that the solver will take it as an identity matrix mathematically. TSAdjoint does not reset the tangent linear solver automatically, TSAdjointResetForward() should be called to reset the tangent linear solver.
955 
956 .seealso: TSSetCostHessianProducts(), TSAdjointResetForward()
957 @*/
958 PetscErrorCode TSAdjointSetForward(TS ts,Mat didp)
959 {
960   Mat            A;
961   Vec            sp;
962   PetscScalar    *xarr;
963   PetscInt       lsize;
964   PetscErrorCode ierr;
965 
966   PetscFunctionBegin;
967   ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */
968   if (!ts->vecs_sensi2) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetCostHessianProducts() first");
969   if (!ts->vec_dir) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Directional vector is missing. Call TSSetCostHessianProducts() to set it.");
970   /* create a single-column dense matrix */
971   ierr = VecGetLocalSize(ts->vec_sol,&lsize);CHKERRQ(ierr);
972   ierr = MatCreateDense(PetscObjectComm((PetscObject)ts),lsize,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&A);CHKERRQ(ierr);
973 
974   ierr = VecDuplicate(ts->vec_sol,&sp);CHKERRQ(ierr);
975   ierr = MatDenseGetColumn(A,0,&xarr);CHKERRQ(ierr);
976   ierr = VecPlaceArray(sp,xarr);CHKERRQ(ierr);
977   if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */
978     if (didp) {
979       ierr = MatMult(didp,ts->vec_dir,sp);CHKERRQ(ierr);
980       ierr = VecScale(sp,2.);CHKERRQ(ierr);
981     } else {
982       ierr = VecZeroEntries(sp);CHKERRQ(ierr);
983     }
984   } else { /* tangent linear variable initialized as dir */
985     ierr = VecCopy(ts->vec_dir,sp);CHKERRQ(ierr);
986   }
987   ierr = VecResetArray(sp);CHKERRQ(ierr);
988   ierr = MatDenseRestoreColumn(A,&xarr);CHKERRQ(ierr);
989   ierr = VecDestroy(&sp);CHKERRQ(ierr);
990 
991   ierr = TSForwardSetInitialSensitivities(ts,A);CHKERRQ(ierr); /* if didp is NULL, identity matrix is assumed */
992 
993   ierr = MatDestroy(&A);CHKERRQ(ierr);
994   PetscFunctionReturn(0);
995 }
996 
997 /*@
998   TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context
999 
1000   Logically Collective on TS
1001 
1002   Input Parameters:
1003 .  ts - the TS context obtained from TSCreate()
1004 
1005   Level: intermediate
1006 
1007 .seealso: TSAdjointSetForward()
1008 @*/
1009 PetscErrorCode TSAdjointResetForward(TS ts)
1010 {
1011   PetscErrorCode ierr;
1012 
1013   PetscFunctionBegin;
1014   ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */
1015   ierr = TSForwardReset(ts);CHKERRQ(ierr);
1016   PetscFunctionReturn(0);
1017 }
1018 
1019 /*@
1020    TSAdjointSetUp - Sets up the internal data structures for the later use
1021    of an adjoint solver
1022 
1023    Collective on TS
1024 
1025    Input Parameter:
1026 .  ts - the TS context obtained from TSCreate()
1027 
1028    Level: advanced
1029 
1030 .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients()
1031 @*/
1032 PetscErrorCode TSAdjointSetUp(TS ts)
1033 {
1034   TSTrajectory     tj;
1035   PetscBool        match;
1036   PetscErrorCode   ierr;
1037 
1038   PetscFunctionBegin;
1039   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1040   if (ts->adjointsetupcalled) PetscFunctionReturn(0);
1041   if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first");
1042   if (ts->vecs_sensip && !ts->Jacp && !ts->Jacprhs) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetRHSJacobianP() or TSSetIJacobianP() first");
1043   ierr = TSGetTrajectory(ts,&tj);CHKERRQ(ierr);
1044   ierr = PetscObjectTypeCompare((PetscObject)tj,TSTRAJECTORYBASIC,&match);CHKERRQ(ierr);
1045   if (match) {
1046     PetscBool solution_only;
1047     ierr = TSTrajectoryGetSolutionOnly(tj,&solution_only);CHKERRQ(ierr);
1048     if (solution_only) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"TSAdjoint cannot use the solution-only mode when choosing the Basic TSTrajectory type. Turn it off with -ts_trajectory_solution_only 0");
1049   }
1050   ierr = TSTrajectorySetUseHistory(tj,PETSC_FALSE);CHKERRQ(ierr); /* not use TSHistory */
1051 
1052   if (ts->quadraturets) { /* if there is integral in the cost function */
1053     ierr = VecDuplicate(ts->vecs_sensi[0],&ts->vec_drdu_col);CHKERRQ(ierr);
1054     if (ts->vecs_sensip){
1055       ierr = VecDuplicate(ts->vecs_sensip[0],&ts->vec_drdp_col);CHKERRQ(ierr);
1056     }
1057   }
1058 
1059   if (ts->ops->adjointsetup) {
1060     ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr);
1061   }
1062   ts->adjointsetupcalled = PETSC_TRUE;
1063   PetscFunctionReturn(0);
1064 }
1065 
1066 /*@
1067    TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats.
1068 
1069    Collective on TS
1070 
1071    Input Parameter:
1072 .  ts - the TS context obtained from TSCreate()
1073 
1074    Level: beginner
1075 
1076 .seealso: TSCreate(), TSAdjointSetUp(), TSADestroy()
1077 @*/
1078 PetscErrorCode TSAdjointReset(TS ts)
1079 {
1080   PetscErrorCode ierr;
1081 
1082   PetscFunctionBegin;
1083   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1084   if (ts->ops->adjointreset) {
1085     ierr = (*ts->ops->adjointreset)(ts);CHKERRQ(ierr);
1086   }
1087   if (ts->quadraturets) { /* if there is integral in the cost function */
1088     ierr = VecDestroy(&ts->vec_drdu_col);CHKERRQ(ierr);
1089     if (ts->vecs_sensip){
1090       ierr = VecDestroy(&ts->vec_drdp_col);CHKERRQ(ierr);
1091     }
1092   }
1093   ts->vecs_sensi         = NULL;
1094   ts->vecs_sensip        = NULL;
1095   ts->vecs_sensi2        = NULL;
1096   ts->vecs_sensi2p       = NULL;
1097   ts->vec_dir            = NULL;
1098   ts->adjointsetupcalled = PETSC_FALSE;
1099   PetscFunctionReturn(0);
1100 }
1101 
1102 /*@
1103    TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time
1104 
1105    Logically Collective on TS
1106 
1107    Input Parameters:
1108 +  ts - the TS context obtained from TSCreate()
1109 -  steps - number of steps to use
1110 
1111    Level: intermediate
1112 
1113    Notes:
1114     Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this
1115           so as to integrate back to less than the original timestep
1116 
1117 .seealso: TSSetExactFinalTime()
1118 @*/
1119 PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps)
1120 {
1121   PetscFunctionBegin;
1122   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1123   PetscValidLogicalCollectiveInt(ts,steps,2);
1124   if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps");
1125   if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps");
1126   ts->adjoint_max_steps = steps;
1127   PetscFunctionReturn(0);
1128 }
1129 
1130 /*@C
1131   TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP()
1132 
1133   Level: deprecated
1134 
1135 @*/
1136 PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
1137 {
1138   PetscErrorCode ierr;
1139 
1140   PetscFunctionBegin;
1141   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1142   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
1143 
1144   ts->rhsjacobianp    = func;
1145   ts->rhsjacobianpctx = ctx;
1146   if(Amat) {
1147     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
1148     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
1149     ts->Jacp = Amat;
1150   }
1151   PetscFunctionReturn(0);
1152 }
1153 
1154 /*@C
1155   TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP()
1156 
1157   Level: deprecated
1158 
1159 @*/
1160 PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat)
1161 {
1162   PetscErrorCode ierr;
1163 
1164   PetscFunctionBegin;
1165   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1166   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
1167   PetscValidPointer(Amat,4);
1168 
1169   PetscStackPush("TS user JacobianP function for sensitivity analysis");
1170   ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr);
1171   PetscStackPop;
1172   PetscFunctionReturn(0);
1173 }
1174 
1175 /*@
1176   TSAdjointComputeDRDYFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian()
1177 
1178   Level: deprecated
1179 
1180 @*/
1181 PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
1182 {
1183   PetscErrorCode ierr;
1184 
1185   PetscFunctionBegin;
1186   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1187   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
1188 
1189   PetscStackPush("TS user DRDY function for sensitivity analysis");
1190   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
1191   PetscStackPop;
1192   PetscFunctionReturn(0);
1193 }
1194 
1195 /*@
1196   TSAdjointComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP()
1197 
1198   Level: deprecated
1199 
1200 @*/
1201 PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
1202 {
1203   PetscErrorCode ierr;
1204 
1205   PetscFunctionBegin;
1206   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1207   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
1208 
1209   PetscStackPush("TS user DRDP function for sensitivity analysis");
1210   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
1211   PetscStackPop;
1212   PetscFunctionReturn(0);
1213 }
1214 
1215 /*@C
1216    TSAdjointMonitorSensi - monitors the first lambda sensitivity
1217 
1218    Level: intermediate
1219 
1220 .seealso: TSAdjointMonitorSet()
1221 @*/
1222 PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
1223 {
1224   PetscErrorCode ierr;
1225   PetscViewer    viewer = vf->viewer;
1226 
1227   PetscFunctionBegin;
1228   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
1229   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
1230   ierr = VecView(lambda[0],viewer);CHKERRQ(ierr);
1231   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1232   PetscFunctionReturn(0);
1233 }
1234 
1235 /*@C
1236    TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
1237 
1238    Collective on TS
1239 
1240    Input Parameters:
1241 +  ts - TS object you wish to monitor
1242 .  name - the monitor type one is seeking
1243 .  help - message indicating what monitoring is done
1244 .  manual - manual page for the monitor
1245 .  monitor - the monitor function
1246 -  monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the TS or PetscViewer objects
1247 
1248    Level: developer
1249 
1250 .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
1251           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
1252           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
1253           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
1254           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
1255           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
1256           PetscOptionsFList(), PetscOptionsEList()
1257 @*/
1258 PetscErrorCode TSAdjointMonitorSetFromOptions(TS ts,const char name[],const char help[], const char manual[],PetscErrorCode (*monitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,PetscViewerAndFormat*),PetscErrorCode (*monitorsetup)(TS,PetscViewerAndFormat*))
1259 {
1260   PetscErrorCode    ierr;
1261   PetscViewer       viewer;
1262   PetscViewerFormat format;
1263   PetscBool         flg;
1264 
1265   PetscFunctionBegin;
1266   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
1267   if (flg) {
1268     PetscViewerAndFormat *vf;
1269     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
1270     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
1271     if (monitorsetup) {
1272       ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr);
1273     }
1274     ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
1275   }
1276   PetscFunctionReturn(0);
1277 }
1278 
1279 /*@C
1280    TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every
1281    timestep to display the iteration's  progress.
1282 
1283    Logically Collective on TS
1284 
1285    Input Parameters:
1286 +  ts - the TS context obtained from TSCreate()
1287 .  adjointmonitor - monitoring routine
1288 .  adjointmctx - [optional] user-defined context for private data for the
1289              monitor routine (use NULL if no context is desired)
1290 -  adjointmonitordestroy - [optional] routine that frees monitor context
1291           (may be NULL)
1292 
1293    Calling sequence of monitor:
1294 $    int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx)
1295 
1296 +    ts - the TS context
1297 .    steps - iteration number (after the final time step the monitor routine is called with a step of -1, this is at the final time which may have
1298                                been interpolated to)
1299 .    time - current time
1300 .    u - current iterate
1301 .    numcost - number of cost functionos
1302 .    lambda - sensitivities to initial conditions
1303 .    mu - sensitivities to parameters
1304 -    adjointmctx - [optional] adjoint monitoring context
1305 
1306    Notes:
1307    This routine adds an additional monitor to the list of monitors that
1308    already has been loaded.
1309 
1310    Fortran Notes:
1311     Only a single monitor function can be set for each TS object
1312 
1313    Level: intermediate
1314 
1315 .seealso: TSAdjointMonitorCancel()
1316 @*/
1317 PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**))
1318 {
1319   PetscErrorCode ierr;
1320   PetscInt       i;
1321   PetscBool      identical;
1322 
1323   PetscFunctionBegin;
1324   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1325   for (i=0; i<ts->numbermonitors;i++) {
1326     ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr);
1327     if (identical) PetscFunctionReturn(0);
1328   }
1329   if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set");
1330   ts->adjointmonitor[ts->numberadjointmonitors]          = adjointmonitor;
1331   ts->adjointmonitordestroy[ts->numberadjointmonitors]   = adjointmdestroy;
1332   ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx;
1333   PetscFunctionReturn(0);
1334 }
1335 
1336 /*@C
1337    TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object.
1338 
1339    Logically Collective on TS
1340 
1341    Input Parameters:
1342 .  ts - the TS context obtained from TSCreate()
1343 
1344    Notes:
1345    There is no way to remove a single, specific monitor.
1346 
1347    Level: intermediate
1348 
1349 .seealso: TSAdjointMonitorSet()
1350 @*/
1351 PetscErrorCode TSAdjointMonitorCancel(TS ts)
1352 {
1353   PetscErrorCode ierr;
1354   PetscInt       i;
1355 
1356   PetscFunctionBegin;
1357   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1358   for (i=0; i<ts->numberadjointmonitors; i++) {
1359     if (ts->adjointmonitordestroy[i]) {
1360       ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
1361     }
1362   }
1363   ts->numberadjointmonitors = 0;
1364   PetscFunctionReturn(0);
1365 }
1366 
1367 /*@C
1368    TSAdjointMonitorDefault - the default monitor of adjoint computations
1369 
1370    Level: intermediate
1371 
1372 .seealso: TSAdjointMonitorSet()
1373 @*/
1374 PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
1375 {
1376   PetscErrorCode ierr;
1377   PetscViewer    viewer = vf->viewer;
1378 
1379   PetscFunctionBegin;
1380   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
1381   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
1382   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1383   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");CHKERRQ(ierr);
1384   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
1385   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
1386   PetscFunctionReturn(0);
1387 }
1388 
1389 /*@C
1390    TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling
1391    VecView() for the sensitivities to initial states at each timestep
1392 
1393    Collective on TS
1394 
1395    Input Parameters:
1396 +  ts - the TS context
1397 .  step - current time-step
1398 .  ptime - current time
1399 .  u - current state
1400 .  numcost - number of cost functions
1401 .  lambda - sensitivities to initial conditions
1402 .  mu - sensitivities to parameters
1403 -  dummy - either a viewer or NULL
1404 
1405    Level: intermediate
1406 
1407 .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView()
1408 @*/
1409 PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy)
1410 {
1411   PetscErrorCode   ierr;
1412   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
1413   PetscDraw        draw;
1414   PetscReal        xl,yl,xr,yr,h;
1415   char             time[32];
1416 
1417   PetscFunctionBegin;
1418   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0);
1419 
1420   ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr);
1421   ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr);
1422   ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr);
1423   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
1424   h    = yl + .95*(yr - yl);
1425   ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr);
1426   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
1427   PetscFunctionReturn(0);
1428 }
1429 
1430 /*
1431    TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options.
1432 
1433    Collective on TSAdjoint
1434 
1435    Input Parameter:
1436 .  ts - the TS context
1437 
1438    Options Database Keys:
1439 +  -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
1440 .  -ts_adjoint_monitor - print information at each adjoint time step
1441 -  -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically
1442 
1443    Level: developer
1444 
1445    Notes:
1446     This is not normally called directly by users
1447 
1448 .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp()
1449 */
1450 PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts)
1451 {
1452   PetscBool      tflg,opt;
1453   PetscErrorCode ierr;
1454 
1455   PetscFunctionBegin;
1456   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
1457   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr);
1458   tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE;
1459   ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr);
1460   if (opt) {
1461     ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
1462     ts->adjoint_solve = tflg;
1463   }
1464   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr);
1465   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr);
1466   opt  = PETSC_FALSE;
1467   ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr);
1468   if (opt) {
1469     TSMonitorDrawCtx ctx;
1470     PetscInt         howoften = 1;
1471 
1472     ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr);
1473     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
1474     ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
1475   }
1476   PetscFunctionReturn(0);
1477 }
1478 
1479 /*@
1480    TSAdjointStep - Steps one time step backward in the adjoint run
1481 
1482    Collective on TS
1483 
1484    Input Parameter:
1485 .  ts - the TS context obtained from TSCreate()
1486 
1487    Level: intermediate
1488 
1489 .seealso: TSAdjointSetUp(), TSAdjointSolve()
1490 @*/
1491 PetscErrorCode TSAdjointStep(TS ts)
1492 {
1493   DM               dm;
1494   PetscErrorCode   ierr;
1495 
1496   PetscFunctionBegin;
1497   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1498   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
1499   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
1500 
1501   ts->reason = TS_CONVERGED_ITERATING;
1502   ts->ptime_prev = ts->ptime;
1503   if (!ts->ops->adjointstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed because the adjoint of  %s has not been implemented, try other time stepping methods for adjoint sensitivity analysis",((PetscObject)ts)->type_name);
1504   ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
1505   ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr);
1506   ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
1507   ts->adjoint_steps++; ts->steps--;
1508 
1509   if (ts->reason < 0) {
1510     if (ts->errorifstepfailed) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]);
1511   } else if (!ts->reason) {
1512     if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1513   }
1514   PetscFunctionReturn(0);
1515 }
1516 
1517 /*@
1518    TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE
1519 
1520    Collective on TS
1521 
1522    Input Parameter:
1523 .  ts - the TS context obtained from TSCreate()
1524 
1525    Options Database:
1526 . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values
1527 
1528    Level: intermediate
1529 
1530    Notes:
1531    This must be called after a call to TSSolve() that solves the forward problem
1532 
1533    By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time
1534 
1535 .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep()
1536 @*/
1537 PetscErrorCode TSAdjointSolve(TS ts)
1538 {
1539 #if defined(TSADJOINT_STAGE)
1540   PetscLogStage  adjoint_stage;
1541 #endif
1542   PetscErrorCode ierr;
1543 
1544   PetscFunctionBegin;
1545   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1546 #if defined(TSADJOINT_STAGE)
1547   ierr = PetscLogStageRegister("TSAdjoint",&adjoint_stage);CHKERRQ(ierr);
1548   ierr = PetscLogStagePush(adjoint_stage);CHKERRQ(ierr);
1549 #endif
1550   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
1551 
1552   /* reset time step and iteration counters */
1553   ts->adjoint_steps     = 0;
1554   ts->ksp_its           = 0;
1555   ts->snes_its          = 0;
1556   ts->num_snes_failures = 0;
1557   ts->reject            = 0;
1558   ts->reason            = TS_CONVERGED_ITERATING;
1559 
1560   if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps;
1561   if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
1562 
1563   while (!ts->reason) {
1564     ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
1565     ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
1566     ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr);
1567     ierr = TSAdjointStep(ts);CHKERRQ(ierr);
1568     if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) {
1569       ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr);
1570     }
1571   }
1572   ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
1573   ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
1574   ts->solvetime = ts->ptime;
1575   ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr);
1576   ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr);
1577   ts->adjoint_max_steps = 0;
1578 #if defined(TSADJOINT_STAGE)
1579   ierr = PetscLogStagePop();CHKERRQ(ierr);
1580 #endif
1581   PetscFunctionReturn(0);
1582 }
1583 
1584 /*@C
1585    TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet()
1586 
1587    Collective on TS
1588 
1589    Input Parameters:
1590 +  ts - time stepping context obtained from TSCreate()
1591 .  step - step number that has just completed
1592 .  ptime - model time of the state
1593 .  u - state at the current model time
1594 .  numcost - number of cost functions (dimension of lambda  or mu)
1595 .  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
1596 -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
1597 
1598    Notes:
1599    TSAdjointMonitor() is typically used automatically within the time stepping implementations.
1600    Users would almost never call this routine directly.
1601 
1602    Level: developer
1603 
1604 @*/
1605 PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu)
1606 {
1607   PetscErrorCode ierr;
1608   PetscInt       i,n = ts->numberadjointmonitors;
1609 
1610   PetscFunctionBegin;
1611   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1612   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
1613   ierr = VecLockReadPush(u);CHKERRQ(ierr);
1614   for (i=0; i<n; i++) {
1615     ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
1616   }
1617   ierr = VecLockReadPop(u);CHKERRQ(ierr);
1618   PetscFunctionReturn(0);
1619 }
1620 
1621 /*@
1622  TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run.
1623 
1624  Collective on TS
1625 
1626  Input Arguments:
1627  .  ts - time stepping context
1628 
1629  Level: advanced
1630 
1631  Notes:
1632  This function cannot be called until TSAdjointStep() has been completed.
1633 
1634  .seealso: TSAdjointSolve(), TSAdjointStep
1635  @*/
1636 PetscErrorCode TSAdjointCostIntegral(TS ts)
1637 {
1638     PetscErrorCode ierr;
1639     PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1640     if (!ts->ops->adjointintegral) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide integral evaluation in the adjoint run",((PetscObject)ts)->type_name);
1641     ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr);
1642     PetscFunctionReturn(0);
1643 }
1644 
1645 /* ------------------ Forward (tangent linear) sensitivity  ------------------*/
1646 
1647 /*@
1648   TSForwardSetUp - Sets up the internal data structures for the later use
1649   of forward sensitivity analysis
1650 
1651   Collective on TS
1652 
1653   Input Parameter:
1654 . ts - the TS context obtained from TSCreate()
1655 
1656   Level: advanced
1657 
1658 .seealso: TSCreate(), TSDestroy(), TSSetUp()
1659 @*/
1660 PetscErrorCode TSForwardSetUp(TS ts)
1661 {
1662   PetscErrorCode ierr;
1663 
1664   PetscFunctionBegin;
1665   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1666   if (ts->forwardsetupcalled) PetscFunctionReturn(0);
1667   if (ts->ops->forwardsetup) {
1668     ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr);
1669   }
1670   ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr);
1671   ts->forwardsetupcalled = PETSC_TRUE;
1672   PetscFunctionReturn(0);
1673 }
1674 
1675 /*@
1676   TSForwardReset - Reset the internal data structures used by forward sensitivity analysis
1677 
1678   Collective on TS
1679 
1680   Input Parameter:
1681 . ts - the TS context obtained from TSCreate()
1682 
1683   Level: advanced
1684 
1685 .seealso: TSCreate(), TSDestroy(), TSForwardSetUp()
1686 @*/
1687 PetscErrorCode TSForwardReset(TS ts)
1688 {
1689   TS             quadts = ts->quadraturets;
1690   PetscErrorCode ierr;
1691 
1692   PetscFunctionBegin;
1693   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1694   if (ts->ops->forwardreset) {
1695     ierr = (*ts->ops->forwardreset)(ts);CHKERRQ(ierr);
1696   }
1697   ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr);
1698   if (quadts) {
1699     ierr = MatDestroy(&quadts->mat_sensip);CHKERRQ(ierr);
1700   }
1701   ierr = VecDestroy(&ts->vec_sensip_col);CHKERRQ(ierr);
1702   ts->forward_solve      = PETSC_FALSE;
1703   ts->forwardsetupcalled = PETSC_FALSE;
1704   PetscFunctionReturn(0);
1705 }
1706 
1707 /*@
1708   TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term.
1709 
1710   Input Parameter:
1711 + ts- the TS context obtained from TSCreate()
1712 . numfwdint- number of integrals
1713 - vp = the vectors containing the gradients for each integral w.r.t. parameters
1714 
1715   Level: deprecated
1716 
1717 .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1718 @*/
1719 PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp)
1720 {
1721   PetscFunctionBegin;
1722   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1723   if (ts->numcost && ts->numcost!=numfwdint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand()");
1724   if (!ts->numcost) ts->numcost = numfwdint;
1725 
1726   ts->vecs_integral_sensip = vp;
1727   PetscFunctionReturn(0);
1728 }
1729 
1730 /*@
1731   TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term.
1732 
1733   Input Parameter:
1734 . ts- the TS context obtained from TSCreate()
1735 
1736   Output Parameter:
1737 . vp = the vectors containing the gradients for each integral w.r.t. parameters
1738 
1739   Level: deprecated
1740 
1741 .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1742 @*/
1743 PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp)
1744 {
1745   PetscFunctionBegin;
1746   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1747   PetscValidPointer(vp,3);
1748   if (numfwdint) *numfwdint = ts->numcost;
1749   if (vp) *vp = ts->vecs_integral_sensip;
1750   PetscFunctionReturn(0);
1751 }
1752 
1753 /*@
1754   TSForwardStep - Compute the forward sensitivity for one time step.
1755 
1756   Collective on TS
1757 
1758   Input Arguments:
1759 . ts - time stepping context
1760 
1761   Level: advanced
1762 
1763   Notes:
1764   This function cannot be called until TSStep() has been completed.
1765 
1766 .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp()
1767 @*/
1768 PetscErrorCode TSForwardStep(TS ts)
1769 {
1770   PetscErrorCode ierr;
1771   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1772   if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name);
1773   ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1774   ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr);
1775   ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1776   if (ts->reason < 0 && ts->errorifstepfailed) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSFowardStep has failed due to %s",TSConvergedReasons[ts->reason]);
1777   PetscFunctionReturn(0);
1778 }
1779 
1780 /*@
1781   TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution  w.r.t. the problem parameters and initial values.
1782 
1783   Logically Collective on TS
1784 
1785   Input Parameters:
1786 + ts - the TS context obtained from TSCreate()
1787 . nump - number of parameters
1788 - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1789 
1790   Level: beginner
1791 
1792   Notes:
1793   Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems.
1794   This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically.
1795   You must call this function before TSSolve().
1796   The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime.
1797 
1798 .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1799 @*/
1800 PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat)
1801 {
1802   PetscErrorCode ierr;
1803 
1804   PetscFunctionBegin;
1805   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1806   PetscValidHeaderSpecific(Smat,MAT_CLASSID,3);
1807   ts->forward_solve  = PETSC_TRUE;
1808   if (nump == PETSC_DEFAULT) {
1809     ierr = MatGetSize(Smat,NULL,&ts->num_parameters);CHKERRQ(ierr);
1810   } else ts->num_parameters = nump;
1811   ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr);
1812   ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr);
1813   ts->mat_sensip = Smat;
1814   PetscFunctionReturn(0);
1815 }
1816 
1817 /*@
1818   TSForwardGetSensitivities - Returns the trajectory sensitivities
1819 
1820   Not Collective, but Vec returned is parallel if TS is parallel
1821 
1822   Output Parameter:
1823 + ts - the TS context obtained from TSCreate()
1824 . nump - number of parameters
1825 - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1826 
1827   Level: intermediate
1828 
1829 .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1830 @*/
1831 PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat)
1832 {
1833   PetscFunctionBegin;
1834   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1835   if (nump) *nump = ts->num_parameters;
1836   if (Smat) *Smat = ts->mat_sensip;
1837   PetscFunctionReturn(0);
1838 }
1839 
1840 /*@
1841    TSForwardCostIntegral - Evaluate the cost integral in the forward run.
1842 
1843    Collective on TS
1844 
1845    Input Arguments:
1846 .  ts - time stepping context
1847 
1848    Level: advanced
1849 
1850    Notes:
1851    This function cannot be called until TSStep() has been completed.
1852 
1853 .seealso: TSSolve(), TSAdjointCostIntegral()
1854 @*/
1855 PetscErrorCode TSForwardCostIntegral(TS ts)
1856 {
1857   PetscErrorCode ierr;
1858   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1859   if (!ts->ops->forwardintegral) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide integral evaluation in the forward run",((PetscObject)ts)->type_name);
1860   ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr);
1861   PetscFunctionReturn(0);
1862 }
1863 
1864 /*@
1865   TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities
1866 
1867   Collective on TS
1868 
1869   Input Parameter
1870 + ts - the TS context obtained from TSCreate()
1871 - didp - parametric sensitivities of the initial condition
1872 
1873   Level: intermediate
1874 
1875   Notes: TSSolve() allows users to pass the initial solution directly to TS. But the tangent linear variables cannot be initialized in this way. This function is used to set initial values for tangent linear variables.
1876 
1877 .seealso: TSForwardSetSensitivities()
1878 @*/
1879 PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp)
1880 {
1881   PetscErrorCode ierr;
1882 
1883   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1884   PetscValidHeaderSpecific(didp,MAT_CLASSID,2);
1885   if (!ts->mat_sensip) {
1886     ierr = TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp);CHKERRQ(ierr);
1887   }
1888   PetscFunctionReturn(0);
1889 }
1890 
1891 /*@
1892    TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages
1893 
1894    Input Parameter:
1895 .  ts - the TS context obtained from TSCreate()
1896 
1897    Output Parameters:
1898 +  ns - number of stages
1899 -  S - tangent linear sensitivities at the intermediate stages
1900 
1901    Level: advanced
1902 
1903 @*/
1904 PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S)
1905 {
1906   PetscErrorCode ierr;
1907 
1908   PetscFunctionBegin;
1909   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
1910 
1911   if (!ts->ops->getstages) *S=NULL;
1912   else {
1913     ierr = (*ts->ops->forwardgetstages)(ts,ns,S);CHKERRQ(ierr);
1914   }
1915   PetscFunctionReturn(0);
1916 }
1917 
1918 /*@
1919    TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time
1920 
1921    Input Parameter:
1922 +  ts - the TS context obtained from TSCreate()
1923 -  fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
1924 
1925    Output Parameters:
1926 .  quadts - the child TS context
1927 
1928    Level: intermediate
1929 
1930 .seealso: TSGetQuadratureTS()
1931 @*/
1932 PetscErrorCode TSCreateQuadratureTS(TS ts,PetscBool fwd,TS *quadts)
1933 {
1934   char prefix[128];
1935   PetscErrorCode ierr;
1936 
1937   PetscFunctionBegin;
1938   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1939   PetscValidPointer(quadts,2);
1940   ierr = TSDestroy(&ts->quadraturets);CHKERRQ(ierr);
1941   ierr = TSCreate(PetscObjectComm((PetscObject)ts),&ts->quadraturets);CHKERRQ(ierr);
1942   ierr = PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets,(PetscObject)ts,1);CHKERRQ(ierr);
1943   ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->quadraturets);CHKERRQ(ierr);
1944   ierr = PetscSNPrintf(prefix,sizeof(prefix),"%squad_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "");CHKERRQ(ierr);
1945   ierr = TSSetOptionsPrefix(ts->quadraturets,prefix);CHKERRQ(ierr);
1946   *quadts = ts->quadraturets;
1947 
1948   if (ts->numcost) {
1949     ierr = VecCreateSeq(PETSC_COMM_SELF,ts->numcost,&(*quadts)->vec_sol);CHKERRQ(ierr);
1950   } else {
1951     ierr = VecCreateSeq(PETSC_COMM_SELF,1,&(*quadts)->vec_sol);CHKERRQ(ierr);
1952   }
1953   ts->costintegralfwd = fwd;
1954   PetscFunctionReturn(0);
1955 }
1956 
1957 /*@
1958    TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time
1959 
1960    Input Parameter:
1961 .  ts - the TS context obtained from TSCreate()
1962 
1963    Output Parameters:
1964 +  fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
1965 -  quadts - the child TS context
1966 
1967    Level: intermediate
1968 
1969 .seealso: TSCreateQuadratureTS()
1970 @*/
1971 PetscErrorCode TSGetQuadratureTS(TS ts,PetscBool *fwd,TS *quadts)
1972 {
1973   PetscFunctionBegin;
1974   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1975   if (fwd) *fwd = ts->costintegralfwd;
1976   if (quadts) *quadts = ts->quadraturets;
1977   PetscFunctionReturn(0);
1978 }
1979 
1980 /*@
1981    TSComputeSNESJacobian - Compute the SNESJacobian
1982 
1983    Input Parameters:
1984 +  ts - the TS context obtained from TSCreate()
1985 -  x - state vector
1986 
1987    Output Parameters:
1988 +  J - Jacobian matrix
1989 -  Jpre - preconditioning matrix for J (may be same as J)
1990 
1991    Level: developer
1992 
1993    Notes:
1994    Using SNES to compute the Jacobian enables finite differencing when TS Jacobian is not available.
1995 @*/
1996 PetscErrorCode TSComputeSNESJacobian(TS ts,Vec x,Mat J,Mat Jpre)
1997 {
1998   SNES           snes = ts->snes;
1999   PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*) = NULL;
2000   PetscErrorCode ierr;
2001 
2002   PetscFunctionBegin;
2003   /*
2004     Unlike implicit methods, explicit methods do not have SNESMatFDColoring in the snes object
2005     because SNESSolve() has not been called yet; so querying SNESMatFDColoring does not work for
2006     explicit methods. Instead, we check the Jacobian compute function directly to determin if FD
2007     coloring is used.
2008   */
2009   ierr = SNESGetJacobian(snes,NULL,NULL,&jac,NULL);CHKERRQ(ierr);
2010   if (jac == SNESComputeJacobianDefaultColor) {
2011     Vec f;
2012     ierr = SNESSetSolution(snes,x);CHKERRQ(ierr);
2013     ierr = SNESGetFunction(snes,&f,NULL,NULL);CHKERRQ(ierr);
2014     /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */
2015     ierr = SNESComputeFunction(snes,x,f);CHKERRQ(ierr);
2016   }
2017   ierr = SNESComputeJacobian(snes,x,J,Jpre);CHKERRQ(ierr);
2018   PetscFunctionReturn(0);
2019 }
2020