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