xref: /petsc/src/ts/interface/sensitivity/tssen.c (revision c9ad9de076a40494a84ae548603ee9f9be9d5015)
1 #include <petsc/private/tsimpl.h>        /*I "petscts.h"  I*/
2 #include <petscdraw.h>
3 
4 PetscLogEvent TS_AdjointStep, TS_ForwardStep;
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   - The TS context obtained from TSCreate()
15 - func - The function
16 
17   Calling sequence of func:
18 $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx);
19 +   t - current timestep
20 .   U - input vector (current ODE solution)
21 .   A - output matrix
22 -   ctx - [optional] user-defined function context
23 
24   Level: intermediate
25 
26   Notes:
27     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
28 
29 .keywords: TS, sensitivity
30 .seealso:
31 @*/
32 PetscErrorCode TSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
33 {
34   PetscErrorCode ierr;
35 
36   PetscFunctionBegin;
37   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
38   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
39 
40   ts->rhsjacobianp    = func;
41   ts->rhsjacobianpctx = ctx;
42   if(Amat) {
43     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
44     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
45     ts->Jacp = Amat;
46   }
47   PetscFunctionReturn(0);
48 }
49 
50 /*@C
51   TSComputeRHSJacobianP - Runs the user-defined JacobianP function.
52 
53   Collective on TS
54 
55   Input Parameters:
56 . ts   - The TS context obtained from TSCreate()
57 
58   Level: developer
59 
60 .keywords: TS, sensitivity
61 .seealso: TSSetRHSJacobianP()
62 @*/
63 PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec U,Mat Amat)
64 {
65   PetscErrorCode ierr;
66 
67   PetscFunctionBegin;
68   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
69   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
70   PetscValidPointer(Amat,4);
71 
72   PetscStackPush("TS user JacobianP function for sensitivity analysis");
73   ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx);CHKERRQ(ierr);
74   PetscStackPop;
75   PetscFunctionReturn(0);
76 }
77 
78 /*@C
79     TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions
80 
81     Logically Collective on TS
82 
83     Input Parameters:
84 +   ts - the TS context obtained from TSCreate()
85 .   numcost - number of gradients to be computed, this is the number of cost functions
86 .   costintegral - vector that stores the integral values
87 .   rf - routine for evaluating the integrand function
88 .   drduf - function that computes the gradients of the r's with respect to u
89 .   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)
90 .   fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run
91 -   ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL)
92 
93     Calling sequence of rf:
94 $   PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx);
95 
96     Calling sequence of drduf:
97 $   PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx);
98 
99     Calling sequence of drdpf:
100 $   PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx);
101 
102     Level: intermediate
103 
104     Notes:
105     For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions
106 
107 .keywords: TS, sensitivity analysis, timestep, set, quadrature, function
108 
109 .seealso: TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients()
110 @*/
111 PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*),
112                                                           PetscErrorCode (*drduf)(TS,PetscReal,Vec,Vec*,void*),
113                                                           PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*),
114                                                           PetscBool fwd,void *ctx)
115 {
116   PetscErrorCode ierr;
117 
118   PetscFunctionBegin;
119   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
120   if (costintegral) PetscValidHeaderSpecific(costintegral,VEC_CLASSID,3);
121   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()");
122   if (!ts->numcost) ts->numcost=numcost;
123 
124   if (costintegral) {
125     ierr = PetscObjectReference((PetscObject)costintegral);CHKERRQ(ierr);
126     ierr = VecDestroy(&ts->vec_costintegral);CHKERRQ(ierr);
127     ts->vec_costintegral = costintegral;
128   } else {
129     if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */
130       ierr = VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);CHKERRQ(ierr);
131     } else {
132       ierr = VecSet(ts->vec_costintegral,0.0);CHKERRQ(ierr);
133     }
134   }
135   if (!ts->vec_costintegrand) {
136     ierr = VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);CHKERRQ(ierr);
137   } else {
138     ierr = VecSet(ts->vec_costintegrand,0.0);CHKERRQ(ierr);
139   }
140   ts->costintegralfwd  = fwd; /* Evaluate the cost integral in forward run if fwd is true */
141   ts->costintegrand    = rf;
142   ts->costintegrandctx = ctx;
143   ts->drdufunction     = drduf;
144   ts->drdpfunction     = drdpf;
145   PetscFunctionReturn(0);
146 }
147 
148 /*@
149    TSGetCostIntegral - Returns the values of the integral term in the cost functions.
150    It is valid to call the routine after a backward run.
151 
152    Not Collective
153 
154    Input Parameter:
155 .  ts - the TS context obtained from TSCreate()
156 
157    Output Parameter:
158 .  v - the vector containing the integrals for each cost function
159 
160    Level: intermediate
161 
162 .seealso: TSSetCostIntegrand()
163 
164 .keywords: TS, sensitivity analysis
165 @*/
166 PetscErrorCode  TSGetCostIntegral(TS ts,Vec *v)
167 {
168   PetscFunctionBegin;
169   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
170   PetscValidPointer(v,2);
171   *v = ts->vec_costintegral;
172   PetscFunctionReturn(0);
173 }
174 
175 /*@
176    TSComputeCostIntegrand - Evaluates the integral function in the cost functions.
177 
178    Input Parameters:
179 +  ts - the TS context
180 .  t - current time
181 -  U - state vector, i.e. current solution
182 
183    Output Parameter:
184 .  Q - vector of size numcost to hold the outputs
185 
186    Note:
187    Most users should not need to explicitly call this routine, as it
188    is used internally within the sensitivity analysis context.
189 
190    Level: developer
191 
192 .keywords: TS, compute
193 
194 .seealso: TSSetCostIntegrand()
195 @*/
196 PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec U,Vec Q)
197 {
198   PetscErrorCode ierr;
199 
200   PetscFunctionBegin;
201   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
202   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
203   PetscValidHeaderSpecific(Q,VEC_CLASSID,4);
204 
205   ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
206   if (ts->costintegrand) {
207     PetscStackPush("TS user integrand in the cost function");
208     ierr = (*ts->costintegrand)(ts,t,U,Q,ts->costintegrandctx);CHKERRQ(ierr);
209     PetscStackPop;
210   } else {
211     ierr = VecZeroEntries(Q);CHKERRQ(ierr);
212   }
213 
214   ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr);
215   PetscFunctionReturn(0);
216 }
217 
218 /*@
219   TSComputeDRDUFunction - Runs the user-defined DRDU function.
220 
221   Collective on TS
222 
223   Input Parameters:
224 + ts - the TS context obtained from TSCreate()
225 . t - current time
226 - U - stata vector
227 
228   Output Parameters:
229 . DRDU - vecotr array to hold the outputs
230 
231   Notes:
232   TSComputeDRDUFunction() is typically used for sensitivity implementation,
233   so most users would not generally call this routine themselves.
234 
235   Level: developer
236 
237 .keywords: TS, sensitivity
238 .seealso: TSSetCostIntegrand()
239 @*/
240 PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
241 {
242   PetscErrorCode ierr;
243 
244   PetscFunctionBegin;
245   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
246   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
247 
248   PetscStackPush("TS user DRDU function for sensitivity analysis");
249   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
250   PetscStackPop;
251   PetscFunctionReturn(0);
252 }
253 
254 /*@
255   TSComputeDRDPFunction - Runs the user-defined DRDP function.
256 
257   Collective on TS
258 
259   Input Parameters:
260 + ts - the TS context obtained from TSCreate()
261 . t - current time
262 - U - stata vector
263 
264   Output Parameters:
265 . DRDP - vecotr array to hold the outputs
266 
267   Notes:
268   TSComputeDRDPFunction() is typically used for sensitivity implementation,
269   so most users would not generally call this routine themselves.
270 
271   Level: developer
272 
273 .keywords: TS, sensitivity
274 .seealso: TSSetCostIntegrand()
275 @*/
276 PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
277 {
278   PetscErrorCode ierr;
279 
280   PetscFunctionBegin;
281   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
282   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
283 
284   PetscStackPush("TS user DRDP function for sensitivity analysis");
285   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
286   PetscStackPop;
287   PetscFunctionReturn(0);
288 }
289 
290 /* --------------------------- Adjoint sensitivity ---------------------------*/
291 
292 /*@
293    TSSetCostGradients - Sets the initial value of the gradients of the cost function w.r.t. initial values and w.r.t. the problem parameters
294       for use by the TSAdjoint routines.
295 
296    Logically Collective on TS and Vec
297 
298    Input Parameters:
299 +  ts - the TS context obtained from TSCreate()
300 .  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
301 -  mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
302 
303    Level: beginner
304 
305    Notes:
306     the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime  mu_i = df/dp|finaltime
307 
308    After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities
309 
310 .keywords: TS, timestep, set, sensitivity, initial values
311 @*/
312 PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu)
313 {
314   PetscFunctionBegin;
315   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
316   PetscValidPointer(lambda,2);
317   ts->vecs_sensi  = lambda;
318   ts->vecs_sensip = mu;
319   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");
320   ts->numcost  = numcost;
321   PetscFunctionReturn(0);
322 }
323 
324 /*@
325    TSGetCostGradients - Returns the gradients from the TSAdjointSolve()
326 
327    Not Collective, but Vec returned is parallel if TS is parallel
328 
329    Input Parameter:
330 .  ts - the TS context obtained from TSCreate()
331 
332    Output Parameter:
333 +  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
334 -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
335 
336    Level: intermediate
337 
338 .seealso: TSGetTimeStep()
339 
340 .keywords: TS, timestep, get, sensitivity
341 @*/
342 PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu)
343 {
344   PetscFunctionBegin;
345   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
346   if (numcost) *numcost = ts->numcost;
347   if (lambda)  *lambda  = ts->vecs_sensi;
348   if (mu)      *mu      = ts->vecs_sensip;
349   PetscFunctionReturn(0);
350 }
351 
352 /*@
353    TSAdjointSetUp - Sets up the internal data structures for the later use
354    of an adjoint solver
355 
356    Collective on TS
357 
358    Input Parameter:
359 .  ts - the TS context obtained from TSCreate()
360 
361    Level: advanced
362 
363 .keywords: TS, timestep, setup
364 
365 .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients()
366 @*/
367 PetscErrorCode TSAdjointSetUp(TS ts)
368 {
369   PetscErrorCode ierr;
370 
371   PetscFunctionBegin;
372   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
373   if (ts->adjointsetupcalled) PetscFunctionReturn(0);
374   if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first");
375   if (ts->vecs_sensip && !ts->Jacp) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSAdjointSetRHSJacobian() first");
376 
377   if (ts->vec_costintegral) { /* if there is integral in the cost function */
378     ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&ts->vecs_drdu);CHKERRQ(ierr);
379     if (ts->vecs_sensip){
380       ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr);
381     }
382   }
383 
384   if (ts->ops->adjointsetup) {
385     ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr);
386   }
387   ts->adjointsetupcalled = PETSC_TRUE;
388   PetscFunctionReturn(0);
389 }
390 
391 /*@
392    TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time
393 
394    Logically Collective on TS
395 
396    Input Parameters:
397 +  ts - the TS context obtained from TSCreate()
398 .  steps - number of steps to use
399 
400    Level: intermediate
401 
402    Notes:
403     Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this
404           so as to integrate back to less than the original timestep
405 
406 .keywords: TS, timestep, set, maximum, iterations
407 
408 .seealso: TSSetExactFinalTime()
409 @*/
410 PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps)
411 {
412   PetscFunctionBegin;
413   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
414   PetscValidLogicalCollectiveInt(ts,steps,2);
415   if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps");
416   if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps");
417   ts->adjoint_max_steps = steps;
418   PetscFunctionReturn(0);
419 }
420 
421 /*@C
422   TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP()
423 
424   Level: deprecated
425 
426 @*/
427 PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx)
428 {
429   PetscErrorCode ierr;
430 
431   PetscFunctionBegin;
432   PetscValidHeaderSpecific(ts, TS_CLASSID,1);
433   PetscValidHeaderSpecific(Amat,MAT_CLASSID,2);
434 
435   ts->rhsjacobianp    = func;
436   ts->rhsjacobianpctx = ctx;
437   if(Amat) {
438     ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr);
439     ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr);
440     ts->Jacp = Amat;
441   }
442   PetscFunctionReturn(0);
443 }
444 
445 /*@C
446   TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP()
447 
448   Level: deprecated
449 
450 @*/
451 PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat)
452 {
453   PetscErrorCode ierr;
454 
455   PetscFunctionBegin;
456   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
457   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
458   PetscValidPointer(Amat,4);
459 
460   PetscStackPush("TS user JacobianP function for sensitivity analysis");
461   ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr);
462   PetscStackPop;
463   PetscFunctionReturn(0);
464 }
465 
466 /*@
467   TSAdjointComputeDRDYFunction - Deprecated, use TSComputeDRDUFunction()
468 
469   Level: deprecated
470 
471 @*/
472 PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU)
473 {
474   PetscErrorCode ierr;
475 
476   PetscFunctionBegin;
477   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
478   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
479 
480   PetscStackPush("TS user DRDY function for sensitivity analysis");
481   ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr);
482   PetscStackPop;
483   PetscFunctionReturn(0);
484 }
485 
486 /*@
487   TSAdjointComputeDRDPFunction - Deprecated, use TSComputeDRDPFunction()
488 
489   Level: deprecated
490 
491 @*/
492 PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP)
493 {
494   PetscErrorCode ierr;
495 
496   PetscFunctionBegin;
497   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
498   PetscValidHeaderSpecific(U,VEC_CLASSID,3);
499 
500   PetscStackPush("TS user DRDP function for sensitivity analysis");
501   ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr);
502   PetscStackPop;
503   PetscFunctionReturn(0);
504 }
505 
506 /*@C
507    TSAdjointMonitorSensi - monitors the first lambda sensitivity
508 
509    Level: intermediate
510 
511 .keywords: TS, set, monitor
512 
513 .seealso: TSAdjointMonitorSet()
514 @*/
515 PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
516 {
517   PetscErrorCode ierr;
518   PetscViewer    viewer = vf->viewer;
519 
520   PetscFunctionBegin;
521   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
522   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
523   ierr = VecView(lambda[0],viewer);CHKERRQ(ierr);
524   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
525   PetscFunctionReturn(0);
526 }
527 
528 /*@C
529    TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user
530 
531    Collective on TS
532 
533    Input Parameters:
534 +  ts - TS object you wish to monitor
535 .  name - the monitor type one is seeking
536 .  help - message indicating what monitoring is done
537 .  manual - manual page for the monitor
538 .  monitor - the monitor function
539 -  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
540 
541    Level: developer
542 
543 .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(),
544           PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool()
545           PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(),
546           PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(),
547           PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(),
548           PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(),
549           PetscOptionsFList(), PetscOptionsEList()
550 @*/
551 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*))
552 {
553   PetscErrorCode    ierr;
554   PetscViewer       viewer;
555   PetscViewerFormat format;
556   PetscBool         flg;
557 
558   PetscFunctionBegin;
559   ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr);
560   if (flg) {
561     PetscViewerAndFormat *vf;
562     ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
563     ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
564     if (monitorsetup) {
565       ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr);
566     }
567     ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr);
568   }
569   PetscFunctionReturn(0);
570 }
571 
572 /*@C
573    TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every
574    timestep to display the iteration's  progress.
575 
576    Logically Collective on TS
577 
578    Input Parameters:
579 +  ts - the TS context obtained from TSCreate()
580 .  adjointmonitor - monitoring routine
581 .  adjointmctx - [optional] user-defined context for private data for the
582              monitor routine (use NULL if no context is desired)
583 -  adjointmonitordestroy - [optional] routine that frees monitor context
584           (may be NULL)
585 
586    Calling sequence of monitor:
587 $    int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx)
588 
589 +    ts - the TS context
590 .    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
591                                been interpolated to)
592 .    time - current time
593 .    u - current iterate
594 .    numcost - number of cost functionos
595 .    lambda - sensitivities to initial conditions
596 .    mu - sensitivities to parameters
597 -    adjointmctx - [optional] adjoint monitoring context
598 
599    Notes:
600    This routine adds an additional monitor to the list of monitors that
601    already has been loaded.
602 
603    Fortran Notes:
604     Only a single monitor function can be set for each TS object
605 
606    Level: intermediate
607 
608 .keywords: TS, timestep, set, adjoint, monitor
609 
610 .seealso: TSAdjointMonitorCancel()
611 @*/
612 PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**))
613 {
614   PetscErrorCode ierr;
615   PetscInt       i;
616   PetscBool      identical;
617 
618   PetscFunctionBegin;
619   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
620   for (i=0; i<ts->numbermonitors;i++) {
621     ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr);
622     if (identical) PetscFunctionReturn(0);
623   }
624   if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set");
625   ts->adjointmonitor[ts->numberadjointmonitors]          = adjointmonitor;
626   ts->adjointmonitordestroy[ts->numberadjointmonitors]   = adjointmdestroy;
627   ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx;
628   PetscFunctionReturn(0);
629 }
630 
631 /*@C
632    TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object.
633 
634    Logically Collective on TS
635 
636    Input Parameters:
637 .  ts - the TS context obtained from TSCreate()
638 
639    Notes:
640    There is no way to remove a single, specific monitor.
641 
642    Level: intermediate
643 
644 .keywords: TS, timestep, set, adjoint, monitor
645 
646 .seealso: TSAdjointMonitorSet()
647 @*/
648 PetscErrorCode TSAdjointMonitorCancel(TS ts)
649 {
650   PetscErrorCode ierr;
651   PetscInt       i;
652 
653   PetscFunctionBegin;
654   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
655   for (i=0; i<ts->numberadjointmonitors; i++) {
656     if (ts->adjointmonitordestroy[i]) {
657       ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
658     }
659   }
660   ts->numberadjointmonitors = 0;
661   PetscFunctionReturn(0);
662 }
663 
664 /*@C
665    TSAdjointMonitorDefault - the default monitor of adjoint computations
666 
667    Level: intermediate
668 
669 .keywords: TS, set, monitor
670 
671 .seealso: TSAdjointMonitorSet()
672 @*/
673 PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf)
674 {
675   PetscErrorCode ierr;
676   PetscViewer    viewer = vf->viewer;
677 
678   PetscFunctionBegin;
679   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4);
680   ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr);
681   ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
682   ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");CHKERRQ(ierr);
683   ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr);
684   ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
685   PetscFunctionReturn(0);
686 }
687 
688 /*@C
689    TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling
690    VecView() for the sensitivities to initial states at each timestep
691 
692    Collective on TS
693 
694    Input Parameters:
695 +  ts - the TS context
696 .  step - current time-step
697 .  ptime - current time
698 .  u - current state
699 .  numcost - number of cost functions
700 .  lambda - sensitivities to initial conditions
701 .  mu - sensitivities to parameters
702 -  dummy - either a viewer or NULL
703 
704    Level: intermediate
705 
706 .keywords: TS,  vector, adjoint, monitor, view
707 
708 .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView()
709 @*/
710 PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy)
711 {
712   PetscErrorCode   ierr;
713   TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy;
714   PetscDraw        draw;
715   PetscReal        xl,yl,xr,yr,h;
716   char             time[32];
717 
718   PetscFunctionBegin;
719   if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0);
720 
721   ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr);
722   ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr);
723   ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr);
724   ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr);
725   h    = yl + .95*(yr - yl);
726   ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr);
727   ierr = PetscDrawFlush(draw);CHKERRQ(ierr);
728   PetscFunctionReturn(0);
729 }
730 
731 /*
732    TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options.
733 
734    Collective on TSAdjoint
735 
736    Input Parameter:
737 .  ts - the TS context
738 
739    Options Database Keys:
740 +  -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory)
741 .  -ts_adjoint_monitor - print information at each adjoint time step
742 -  -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically
743 
744    Level: developer
745 
746    Notes:
747     This is not normally called directly by users
748 
749 .keywords: TS, trajectory, timestep, set, options, database
750 
751 .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp()
752 */
753 PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts)
754 {
755   PetscBool      tflg,opt;
756   PetscErrorCode ierr;
757 
758   PetscFunctionBegin;
759   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
760   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr);
761   tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE;
762   ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr);
763   if (opt) {
764     ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
765     ts->adjoint_solve = tflg;
766   }
767   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr);
768   ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr);
769   opt  = PETSC_FALSE;
770   ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr);
771   if (opt) {
772     TSMonitorDrawCtx ctx;
773     PetscInt         howoften = 1;
774 
775     ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr);
776     ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr);
777     ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr);
778   }
779   PetscFunctionReturn(0);
780 }
781 
782 /*@
783    TSAdjointStep - Steps one time step backward in the adjoint run
784 
785    Collective on TS
786 
787    Input Parameter:
788 .  ts - the TS context obtained from TSCreate()
789 
790    Level: intermediate
791 
792 .keywords: TS, adjoint, step
793 
794 .seealso: TSAdjointSetUp(), TSAdjointSolve()
795 @*/
796 PetscErrorCode TSAdjointStep(TS ts)
797 {
798   DM               dm;
799   PetscErrorCode   ierr;
800 
801   PetscFunctionBegin;
802   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
803   ierr = TSGetDM(ts,&dm);CHKERRQ(ierr);
804   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
805 
806   ierr = VecViewFromOptions(ts->vec_sol,(PetscObject)ts,"-ts_view_solution");CHKERRQ(ierr);
807 
808   ts->reason = TS_CONVERGED_ITERATING;
809   ts->ptime_prev = ts->ptime;
810   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);
811   ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
812   ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr);
813   ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr);
814   ts->adjoint_steps++; ts->steps--;
815 
816   if (ts->reason < 0) {
817     if (ts->errorifstepfailed) {
818       if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_snes_failures or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
819       else if (ts->reason == TS_DIVERGED_STEP_REJECTED) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_reject or make negative to attempt recovery",TSConvergedReasons[ts->reason]);
820       else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]);
821     }
822   } else if (!ts->reason) {
823     if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
824   }
825   PetscFunctionReturn(0);
826 }
827 
828 /*@
829    TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE
830 
831    Collective on TS
832 
833    Input Parameter:
834 .  ts - the TS context obtained from TSCreate()
835 
836    Options Database:
837 . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values
838 
839    Level: intermediate
840 
841    Notes:
842    This must be called after a call to TSSolve() that solves the forward problem
843 
844    By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time
845 
846 .keywords: TS, timestep, solve
847 
848 .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep()
849 @*/
850 PetscErrorCode TSAdjointSolve(TS ts)
851 {
852   PetscErrorCode    ierr;
853 
854   PetscFunctionBegin;
855   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
856   ierr = TSAdjointSetUp(ts);CHKERRQ(ierr);
857 
858   /* reset time step and iteration counters */
859   ts->adjoint_steps     = 0;
860   ts->ksp_its           = 0;
861   ts->snes_its          = 0;
862   ts->num_snes_failures = 0;
863   ts->reject            = 0;
864   ts->reason            = TS_CONVERGED_ITERATING;
865 
866   if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps;
867   if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS;
868 
869   while (!ts->reason) {
870     ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
871     ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
872     ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr);
873     ierr = TSAdjointStep(ts);CHKERRQ(ierr);
874     if (ts->vec_costintegral && !ts->costintegralfwd) {
875       ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr);
876     }
877   }
878   ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr);
879   ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr);
880   ts->solvetime = ts->ptime;
881   ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr);
882   ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr);
883   ts->adjoint_max_steps = 0;
884   PetscFunctionReturn(0);
885 }
886 
887 /*@C
888    TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet()
889 
890    Collective on TS
891 
892    Input Parameters:
893 +  ts - time stepping context obtained from TSCreate()
894 .  step - step number that has just completed
895 .  ptime - model time of the state
896 .  u - state at the current model time
897 .  numcost - number of cost functions (dimension of lambda  or mu)
898 .  lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables
899 -  mu - vectors containing the gradients of the cost functions with respect to the problem parameters
900 
901    Notes:
902    TSAdjointMonitor() is typically used automatically within the time stepping implementations.
903    Users would almost never call this routine directly.
904 
905    Level: developer
906 
907 .keywords: TS, timestep
908 @*/
909 PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu)
910 {
911   PetscErrorCode ierr;
912   PetscInt       i,n = ts->numberadjointmonitors;
913 
914   PetscFunctionBegin;
915   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
916   PetscValidHeaderSpecific(u,VEC_CLASSID,4);
917   ierr = VecLockReadPush(u);CHKERRQ(ierr);
918   for (i=0; i<n; i++) {
919     ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr);
920   }
921   ierr = VecLockReadPop(u);CHKERRQ(ierr);
922   PetscFunctionReturn(0);
923 }
924 
925 /*@
926  TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run.
927 
928  Collective on TS
929 
930  Input Arguments:
931  .  ts - time stepping context
932 
933  Level: advanced
934 
935  Notes:
936  This function cannot be called until TSAdjointStep() has been completed.
937 
938  .seealso: TSAdjointSolve(), TSAdjointStep
939  @*/
940 PetscErrorCode TSAdjointCostIntegral(TS ts)
941 {
942     PetscErrorCode ierr;
943     PetscValidHeaderSpecific(ts,TS_CLASSID,1);
944     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);
945     ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr);
946     PetscFunctionReturn(0);
947 }
948 
949 /* ------------------ Forward (tangent linear) sensitivity  ------------------*/
950 
951 /*@
952   TSForwardSetUp - Sets up the internal data structures for the later use
953   of forward sensitivity analysis
954 
955   Collective on TS
956 
957   Input Parameter:
958 . ts - the TS context obtained from TSCreate()
959 
960   Level: advanced
961 
962 .keywords: TS, forward sensitivity, setup
963 
964 .seealso: TSCreate(), TSDestroy(), TSSetUp()
965 @*/
966 PetscErrorCode TSForwardSetUp(TS ts)
967 {
968   PetscErrorCode ierr;
969 
970   PetscFunctionBegin;
971   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
972   if (ts->forwardsetupcalled) PetscFunctionReturn(0);
973   if (ts->vec_costintegral && !ts->vecs_integral_sensip ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSForwardSetIntegralGradients() before TSSetCostIntegrand()");
974   if (ts->vecs_integral_sensip) {
975     ierr = VecDuplicateVecs(ts->vec_sol,ts->numcost,&ts->vecs_drdu);CHKERRQ(ierr);
976     ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr);
977   }
978 
979   if (ts->ops->forwardsetup) {
980     ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr);
981   }
982   ts->forwardsetupcalled = PETSC_TRUE;
983   PetscFunctionReturn(0);
984 }
985 
986 /*@
987   TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term.
988 
989   Input Parameter:
990 . ts- the TS context obtained from TSCreate()
991 . numfwdint- number of integrals
992 . vp = the vectors containing the gradients for each integral w.r.t. parameters
993 
994   Level: intermediate
995 
996 .keywords: TS, forward sensitivity
997 
998 .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
999 @*/
1000 PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp)
1001 {
1002   PetscFunctionBegin;
1003   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1004   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()");
1005   if (!ts->numcost) ts->numcost = numfwdint;
1006 
1007   ts->vecs_integral_sensip = vp;
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 /*@
1012   TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term.
1013 
1014   Input Parameter:
1015 . ts- the TS context obtained from TSCreate()
1016 
1017   Output Parameter:
1018 . vp = the vectors containing the gradients for each integral w.r.t. parameters
1019 
1020   Level: intermediate
1021 
1022 .keywords: TS, forward sensitivity
1023 
1024 .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1025 @*/
1026 PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp)
1027 {
1028   PetscFunctionBegin;
1029   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1030   PetscValidPointer(vp,3);
1031   if (numfwdint) *numfwdint = ts->numcost;
1032   if (vp) *vp = ts->vecs_integral_sensip;
1033   PetscFunctionReturn(0);
1034 }
1035 
1036 /*@
1037   TSForwardStep - Compute the forward sensitivity for one time step.
1038 
1039   Collective on TS
1040 
1041   Input Arguments:
1042 . ts - time stepping context
1043 
1044   Level: advanced
1045 
1046   Notes:
1047   This function cannot be called until TSStep() has been completed.
1048 
1049 .keywords: TS, forward sensitivity
1050 
1051 .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp()
1052 @*/
1053 PetscErrorCode TSForwardStep(TS ts)
1054 {
1055   PetscErrorCode ierr;
1056   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1057   if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name);
1058   ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1059   ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr);
1060   ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr);
1061   PetscFunctionReturn(0);
1062 }
1063 
1064 /*@
1065   TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution  w.r.t. the problem parameters and initial values.
1066 
1067   Logically Collective on TS and Vec
1068 
1069   Input Parameters:
1070 + ts - the TS context obtained from TSCreate()
1071 . nump - number of parameters
1072 - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1073 
1074   Level: beginner
1075 
1076   Notes:
1077   Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems.
1078   This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically.
1079   You must call this function before TSSolve().
1080   The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime.
1081 
1082 .keywords: TS, timestep, set, forward sensitivity, initial values
1083 
1084 .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1085 @*/
1086 PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat)
1087 {
1088   PetscErrorCode ierr;
1089 
1090   PetscFunctionBegin;
1091   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1092   PetscValidHeaderSpecific(Smat,MAT_CLASSID,3);
1093   ts->forward_solve  = PETSC_TRUE;
1094   ts->num_parameters = nump;
1095   ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr);
1096   ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr);
1097   ts->mat_sensip = Smat;
1098   PetscFunctionReturn(0);
1099 }
1100 
1101 /*@
1102   TSForwardGetSensitivities - Returns the trajectory sensitivities
1103 
1104   Not Collective, but Vec returned is parallel if TS is parallel
1105 
1106   Output Parameter:
1107 + ts - the TS context obtained from TSCreate()
1108 . nump - number of parameters
1109 - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters
1110 
1111   Level: intermediate
1112 
1113 .keywords: TS, forward sensitivity
1114 
1115 .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep()
1116 @*/
1117 PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat)
1118 {
1119   PetscFunctionBegin;
1120   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1121   if (nump) *nump = ts->num_parameters;
1122   if (Smat) *Smat = ts->mat_sensip;
1123   PetscFunctionReturn(0);
1124 }
1125 
1126 /*@
1127    TSForwardCostIntegral - Evaluate the cost integral in the forward run.
1128 
1129    Collective on TS
1130 
1131    Input Arguments:
1132 .  ts - time stepping context
1133 
1134    Level: advanced
1135 
1136    Notes:
1137    This function cannot be called until TSStep() has been completed.
1138 
1139 .seealso: TSSolve(), TSAdjointCostIntegral()
1140 @*/
1141 PetscErrorCode TSForwardCostIntegral(TS ts)
1142 {
1143   PetscErrorCode ierr;
1144   PetscValidHeaderSpecific(ts,TS_CLASSID,1);
1145   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);
1146   ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr);
1147   PetscFunctionReturn(0);
1148 }
1149