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