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