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